Chris@87
|
1 """
|
Chris@87
|
2 Objects for dealing with polynomials.
|
Chris@87
|
3
|
Chris@87
|
4 This module provides a number of objects (mostly functions) useful for
|
Chris@87
|
5 dealing with polynomials, including a `Polynomial` class that
|
Chris@87
|
6 encapsulates the usual arithmetic operations. (General information
|
Chris@87
|
7 on how this module represents and works with polynomial objects is in
|
Chris@87
|
8 the docstring for its "parent" sub-package, `numpy.polynomial`).
|
Chris@87
|
9
|
Chris@87
|
10 Constants
|
Chris@87
|
11 ---------
|
Chris@87
|
12 - `polydomain` -- Polynomial default domain, [-1,1].
|
Chris@87
|
13 - `polyzero` -- (Coefficients of the) "zero polynomial."
|
Chris@87
|
14 - `polyone` -- (Coefficients of the) constant polynomial 1.
|
Chris@87
|
15 - `polyx` -- (Coefficients of the) identity map polynomial, ``f(x) = x``.
|
Chris@87
|
16
|
Chris@87
|
17 Arithmetic
|
Chris@87
|
18 ----------
|
Chris@87
|
19 - `polyadd` -- add two polynomials.
|
Chris@87
|
20 - `polysub` -- subtract one polynomial from another.
|
Chris@87
|
21 - `polymul` -- multiply two polynomials.
|
Chris@87
|
22 - `polydiv` -- divide one polynomial by another.
|
Chris@87
|
23 - `polypow` -- raise a polynomial to an positive integer power
|
Chris@87
|
24 - `polyval` -- evaluate a polynomial at given points.
|
Chris@87
|
25 - `polyval2d` -- evaluate a 2D polynomial at given points.
|
Chris@87
|
26 - `polyval3d` -- evaluate a 3D polynomial at given points.
|
Chris@87
|
27 - `polygrid2d` -- evaluate a 2D polynomial on a Cartesian product.
|
Chris@87
|
28 - `polygrid3d` -- evaluate a 3D polynomial on a Cartesian product.
|
Chris@87
|
29
|
Chris@87
|
30 Calculus
|
Chris@87
|
31 --------
|
Chris@87
|
32 - `polyder` -- differentiate a polynomial.
|
Chris@87
|
33 - `polyint` -- integrate a polynomial.
|
Chris@87
|
34
|
Chris@87
|
35 Misc Functions
|
Chris@87
|
36 --------------
|
Chris@87
|
37 - `polyfromroots` -- create a polynomial with specified roots.
|
Chris@87
|
38 - `polyroots` -- find the roots of a polynomial.
|
Chris@87
|
39 - `polyvander` -- Vandermonde-like matrix for powers.
|
Chris@87
|
40 - `polyvander2d` -- Vandermonde-like matrix for 2D power series.
|
Chris@87
|
41 - `polyvander3d` -- Vandermonde-like matrix for 3D power series.
|
Chris@87
|
42 - `polycompanion` -- companion matrix in power series form.
|
Chris@87
|
43 - `polyfit` -- least-squares fit returning a polynomial.
|
Chris@87
|
44 - `polytrim` -- trim leading coefficients from a polynomial.
|
Chris@87
|
45 - `polyline` -- polynomial representing given straight line.
|
Chris@87
|
46
|
Chris@87
|
47 Classes
|
Chris@87
|
48 -------
|
Chris@87
|
49 - `Polynomial` -- polynomial class.
|
Chris@87
|
50
|
Chris@87
|
51 See Also
|
Chris@87
|
52 --------
|
Chris@87
|
53 `numpy.polynomial`
|
Chris@87
|
54
|
Chris@87
|
55 """
|
Chris@87
|
56 from __future__ import division, absolute_import, print_function
|
Chris@87
|
57
|
Chris@87
|
58 __all__ = [
|
Chris@87
|
59 'polyzero', 'polyone', 'polyx', 'polydomain', 'polyline', 'polyadd',
|
Chris@87
|
60 'polysub', 'polymulx', 'polymul', 'polydiv', 'polypow', 'polyval',
|
Chris@87
|
61 'polyder', 'polyint', 'polyfromroots', 'polyvander', 'polyfit',
|
Chris@87
|
62 'polytrim', 'polyroots', 'Polynomial', 'polyval2d', 'polyval3d',
|
Chris@87
|
63 'polygrid2d', 'polygrid3d', 'polyvander2d', 'polyvander3d']
|
Chris@87
|
64
|
Chris@87
|
65 import warnings
|
Chris@87
|
66 import numpy as np
|
Chris@87
|
67 import numpy.linalg as la
|
Chris@87
|
68
|
Chris@87
|
69 from . import polyutils as pu
|
Chris@87
|
70 from ._polybase import ABCPolyBase
|
Chris@87
|
71
|
Chris@87
|
72 polytrim = pu.trimcoef
|
Chris@87
|
73
|
Chris@87
|
74 #
|
Chris@87
|
75 # These are constant arrays are of integer type so as to be compatible
|
Chris@87
|
76 # with the widest range of other types, such as Decimal.
|
Chris@87
|
77 #
|
Chris@87
|
78
|
Chris@87
|
79 # Polynomial default domain.
|
Chris@87
|
80 polydomain = np.array([-1, 1])
|
Chris@87
|
81
|
Chris@87
|
82 # Polynomial coefficients representing zero.
|
Chris@87
|
83 polyzero = np.array([0])
|
Chris@87
|
84
|
Chris@87
|
85 # Polynomial coefficients representing one.
|
Chris@87
|
86 polyone = np.array([1])
|
Chris@87
|
87
|
Chris@87
|
88 # Polynomial coefficients representing the identity x.
|
Chris@87
|
89 polyx = np.array([0, 1])
|
Chris@87
|
90
|
Chris@87
|
91 #
|
Chris@87
|
92 # Polynomial series functions
|
Chris@87
|
93 #
|
Chris@87
|
94
|
Chris@87
|
95
|
Chris@87
|
96 def polyline(off, scl):
|
Chris@87
|
97 """
|
Chris@87
|
98 Returns an array representing a linear polynomial.
|
Chris@87
|
99
|
Chris@87
|
100 Parameters
|
Chris@87
|
101 ----------
|
Chris@87
|
102 off, scl : scalars
|
Chris@87
|
103 The "y-intercept" and "slope" of the line, respectively.
|
Chris@87
|
104
|
Chris@87
|
105 Returns
|
Chris@87
|
106 -------
|
Chris@87
|
107 y : ndarray
|
Chris@87
|
108 This module's representation of the linear polynomial ``off +
|
Chris@87
|
109 scl*x``.
|
Chris@87
|
110
|
Chris@87
|
111 See Also
|
Chris@87
|
112 --------
|
Chris@87
|
113 chebline
|
Chris@87
|
114
|
Chris@87
|
115 Examples
|
Chris@87
|
116 --------
|
Chris@87
|
117 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
118 >>> P.polyline(1,-1)
|
Chris@87
|
119 array([ 1, -1])
|
Chris@87
|
120 >>> P.polyval(1, P.polyline(1,-1)) # should be 0
|
Chris@87
|
121 0.0
|
Chris@87
|
122
|
Chris@87
|
123 """
|
Chris@87
|
124 if scl != 0:
|
Chris@87
|
125 return np.array([off, scl])
|
Chris@87
|
126 else:
|
Chris@87
|
127 return np.array([off])
|
Chris@87
|
128
|
Chris@87
|
129
|
Chris@87
|
130 def polyfromroots(roots):
|
Chris@87
|
131 """
|
Chris@87
|
132 Generate a monic polynomial with given roots.
|
Chris@87
|
133
|
Chris@87
|
134 Return the coefficients of the polynomial
|
Chris@87
|
135
|
Chris@87
|
136 .. math:: p(x) = (x - r_0) * (x - r_1) * ... * (x - r_n),
|
Chris@87
|
137
|
Chris@87
|
138 where the `r_n` are the roots specified in `roots`. If a zero has
|
Chris@87
|
139 multiplicity n, then it must appear in `roots` n times. For instance,
|
Chris@87
|
140 if 2 is a root of multiplicity three and 3 is a root of multiplicity 2,
|
Chris@87
|
141 then `roots` looks something like [2, 2, 2, 3, 3]. The roots can appear
|
Chris@87
|
142 in any order.
|
Chris@87
|
143
|
Chris@87
|
144 If the returned coefficients are `c`, then
|
Chris@87
|
145
|
Chris@87
|
146 .. math:: p(x) = c_0 + c_1 * x + ... + x^n
|
Chris@87
|
147
|
Chris@87
|
148 The coefficient of the last term is 1 for monic polynomials in this
|
Chris@87
|
149 form.
|
Chris@87
|
150
|
Chris@87
|
151 Parameters
|
Chris@87
|
152 ----------
|
Chris@87
|
153 roots : array_like
|
Chris@87
|
154 Sequence containing the roots.
|
Chris@87
|
155
|
Chris@87
|
156 Returns
|
Chris@87
|
157 -------
|
Chris@87
|
158 out : ndarray
|
Chris@87
|
159 1-D array of the polynomial's coefficients If all the roots are
|
Chris@87
|
160 real, then `out` is also real, otherwise it is complex. (see
|
Chris@87
|
161 Examples below).
|
Chris@87
|
162
|
Chris@87
|
163 See Also
|
Chris@87
|
164 --------
|
Chris@87
|
165 chebfromroots, legfromroots, lagfromroots, hermfromroots
|
Chris@87
|
166 hermefromroots
|
Chris@87
|
167
|
Chris@87
|
168 Notes
|
Chris@87
|
169 -----
|
Chris@87
|
170 The coefficients are determined by multiplying together linear factors
|
Chris@87
|
171 of the form `(x - r_i)`, i.e.
|
Chris@87
|
172
|
Chris@87
|
173 .. math:: p(x) = (x - r_0) (x - r_1) ... (x - r_n)
|
Chris@87
|
174
|
Chris@87
|
175 where ``n == len(roots) - 1``; note that this implies that `1` is always
|
Chris@87
|
176 returned for :math:`a_n`.
|
Chris@87
|
177
|
Chris@87
|
178 Examples
|
Chris@87
|
179 --------
|
Chris@87
|
180 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
181 >>> P.polyfromroots((-1,0,1)) # x(x - 1)(x + 1) = x^3 - x
|
Chris@87
|
182 array([ 0., -1., 0., 1.])
|
Chris@87
|
183 >>> j = complex(0,1)
|
Chris@87
|
184 >>> P.polyfromroots((-j,j)) # complex returned, though values are real
|
Chris@87
|
185 array([ 1.+0.j, 0.+0.j, 1.+0.j])
|
Chris@87
|
186
|
Chris@87
|
187 """
|
Chris@87
|
188 if len(roots) == 0:
|
Chris@87
|
189 return np.ones(1)
|
Chris@87
|
190 else:
|
Chris@87
|
191 [roots] = pu.as_series([roots], trim=False)
|
Chris@87
|
192 roots.sort()
|
Chris@87
|
193 p = [polyline(-r, 1) for r in roots]
|
Chris@87
|
194 n = len(p)
|
Chris@87
|
195 while n > 1:
|
Chris@87
|
196 m, r = divmod(n, 2)
|
Chris@87
|
197 tmp = [polymul(p[i], p[i+m]) for i in range(m)]
|
Chris@87
|
198 if r:
|
Chris@87
|
199 tmp[0] = polymul(tmp[0], p[-1])
|
Chris@87
|
200 p = tmp
|
Chris@87
|
201 n = m
|
Chris@87
|
202 return p[0]
|
Chris@87
|
203
|
Chris@87
|
204
|
Chris@87
|
205 def polyadd(c1, c2):
|
Chris@87
|
206 """
|
Chris@87
|
207 Add one polynomial to another.
|
Chris@87
|
208
|
Chris@87
|
209 Returns the sum of two polynomials `c1` + `c2`. The arguments are
|
Chris@87
|
210 sequences of coefficients from lowest order term to highest, i.e.,
|
Chris@87
|
211 [1,2,3] represents the polynomial ``1 + 2*x + 3*x**2``.
|
Chris@87
|
212
|
Chris@87
|
213 Parameters
|
Chris@87
|
214 ----------
|
Chris@87
|
215 c1, c2 : array_like
|
Chris@87
|
216 1-D arrays of polynomial coefficients ordered from low to high.
|
Chris@87
|
217
|
Chris@87
|
218 Returns
|
Chris@87
|
219 -------
|
Chris@87
|
220 out : ndarray
|
Chris@87
|
221 The coefficient array representing their sum.
|
Chris@87
|
222
|
Chris@87
|
223 See Also
|
Chris@87
|
224 --------
|
Chris@87
|
225 polysub, polymul, polydiv, polypow
|
Chris@87
|
226
|
Chris@87
|
227 Examples
|
Chris@87
|
228 --------
|
Chris@87
|
229 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
230 >>> c1 = (1,2,3)
|
Chris@87
|
231 >>> c2 = (3,2,1)
|
Chris@87
|
232 >>> sum = P.polyadd(c1,c2); sum
|
Chris@87
|
233 array([ 4., 4., 4.])
|
Chris@87
|
234 >>> P.polyval(2, sum) # 4 + 4(2) + 4(2**2)
|
Chris@87
|
235 28.0
|
Chris@87
|
236
|
Chris@87
|
237 """
|
Chris@87
|
238 # c1, c2 are trimmed copies
|
Chris@87
|
239 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
240 if len(c1) > len(c2):
|
Chris@87
|
241 c1[:c2.size] += c2
|
Chris@87
|
242 ret = c1
|
Chris@87
|
243 else:
|
Chris@87
|
244 c2[:c1.size] += c1
|
Chris@87
|
245 ret = c2
|
Chris@87
|
246 return pu.trimseq(ret)
|
Chris@87
|
247
|
Chris@87
|
248
|
Chris@87
|
249 def polysub(c1, c2):
|
Chris@87
|
250 """
|
Chris@87
|
251 Subtract one polynomial from another.
|
Chris@87
|
252
|
Chris@87
|
253 Returns the difference of two polynomials `c1` - `c2`. The arguments
|
Chris@87
|
254 are sequences of coefficients from lowest order term to highest, i.e.,
|
Chris@87
|
255 [1,2,3] represents the polynomial ``1 + 2*x + 3*x**2``.
|
Chris@87
|
256
|
Chris@87
|
257 Parameters
|
Chris@87
|
258 ----------
|
Chris@87
|
259 c1, c2 : array_like
|
Chris@87
|
260 1-D arrays of polynomial coefficients ordered from low to
|
Chris@87
|
261 high.
|
Chris@87
|
262
|
Chris@87
|
263 Returns
|
Chris@87
|
264 -------
|
Chris@87
|
265 out : ndarray
|
Chris@87
|
266 Of coefficients representing their difference.
|
Chris@87
|
267
|
Chris@87
|
268 See Also
|
Chris@87
|
269 --------
|
Chris@87
|
270 polyadd, polymul, polydiv, polypow
|
Chris@87
|
271
|
Chris@87
|
272 Examples
|
Chris@87
|
273 --------
|
Chris@87
|
274 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
275 >>> c1 = (1,2,3)
|
Chris@87
|
276 >>> c2 = (3,2,1)
|
Chris@87
|
277 >>> P.polysub(c1,c2)
|
Chris@87
|
278 array([-2., 0., 2.])
|
Chris@87
|
279 >>> P.polysub(c2,c1) # -P.polysub(c1,c2)
|
Chris@87
|
280 array([ 2., 0., -2.])
|
Chris@87
|
281
|
Chris@87
|
282 """
|
Chris@87
|
283 # c1, c2 are trimmed copies
|
Chris@87
|
284 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
285 if len(c1) > len(c2):
|
Chris@87
|
286 c1[:c2.size] -= c2
|
Chris@87
|
287 ret = c1
|
Chris@87
|
288 else:
|
Chris@87
|
289 c2 = -c2
|
Chris@87
|
290 c2[:c1.size] += c1
|
Chris@87
|
291 ret = c2
|
Chris@87
|
292 return pu.trimseq(ret)
|
Chris@87
|
293
|
Chris@87
|
294
|
Chris@87
|
295 def polymulx(c):
|
Chris@87
|
296 """Multiply a polynomial by x.
|
Chris@87
|
297
|
Chris@87
|
298 Multiply the polynomial `c` by x, where x is the independent
|
Chris@87
|
299 variable.
|
Chris@87
|
300
|
Chris@87
|
301
|
Chris@87
|
302 Parameters
|
Chris@87
|
303 ----------
|
Chris@87
|
304 c : array_like
|
Chris@87
|
305 1-D array of polynomial coefficients ordered from low to
|
Chris@87
|
306 high.
|
Chris@87
|
307
|
Chris@87
|
308 Returns
|
Chris@87
|
309 -------
|
Chris@87
|
310 out : ndarray
|
Chris@87
|
311 Array representing the result of the multiplication.
|
Chris@87
|
312
|
Chris@87
|
313 Notes
|
Chris@87
|
314 -----
|
Chris@87
|
315
|
Chris@87
|
316 .. versionadded:: 1.5.0
|
Chris@87
|
317
|
Chris@87
|
318 """
|
Chris@87
|
319 # c is a trimmed copy
|
Chris@87
|
320 [c] = pu.as_series([c])
|
Chris@87
|
321 # The zero series needs special treatment
|
Chris@87
|
322 if len(c) == 1 and c[0] == 0:
|
Chris@87
|
323 return c
|
Chris@87
|
324
|
Chris@87
|
325 prd = np.empty(len(c) + 1, dtype=c.dtype)
|
Chris@87
|
326 prd[0] = c[0]*0
|
Chris@87
|
327 prd[1:] = c
|
Chris@87
|
328 return prd
|
Chris@87
|
329
|
Chris@87
|
330
|
Chris@87
|
331 def polymul(c1, c2):
|
Chris@87
|
332 """
|
Chris@87
|
333 Multiply one polynomial by another.
|
Chris@87
|
334
|
Chris@87
|
335 Returns the product of two polynomials `c1` * `c2`. The arguments are
|
Chris@87
|
336 sequences of coefficients, from lowest order term to highest, e.g.,
|
Chris@87
|
337 [1,2,3] represents the polynomial ``1 + 2*x + 3*x**2.``
|
Chris@87
|
338
|
Chris@87
|
339 Parameters
|
Chris@87
|
340 ----------
|
Chris@87
|
341 c1, c2 : array_like
|
Chris@87
|
342 1-D arrays of coefficients representing a polynomial, relative to the
|
Chris@87
|
343 "standard" basis, and ordered from lowest order term to highest.
|
Chris@87
|
344
|
Chris@87
|
345 Returns
|
Chris@87
|
346 -------
|
Chris@87
|
347 out : ndarray
|
Chris@87
|
348 Of the coefficients of their product.
|
Chris@87
|
349
|
Chris@87
|
350 See Also
|
Chris@87
|
351 --------
|
Chris@87
|
352 polyadd, polysub, polydiv, polypow
|
Chris@87
|
353
|
Chris@87
|
354 Examples
|
Chris@87
|
355 --------
|
Chris@87
|
356 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
357 >>> c1 = (1,2,3)
|
Chris@87
|
358 >>> c2 = (3,2,1)
|
Chris@87
|
359 >>> P.polymul(c1,c2)
|
Chris@87
|
360 array([ 3., 8., 14., 8., 3.])
|
Chris@87
|
361
|
Chris@87
|
362 """
|
Chris@87
|
363 # c1, c2 are trimmed copies
|
Chris@87
|
364 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
365 ret = np.convolve(c1, c2)
|
Chris@87
|
366 return pu.trimseq(ret)
|
Chris@87
|
367
|
Chris@87
|
368
|
Chris@87
|
369 def polydiv(c1, c2):
|
Chris@87
|
370 """
|
Chris@87
|
371 Divide one polynomial by another.
|
Chris@87
|
372
|
Chris@87
|
373 Returns the quotient-with-remainder of two polynomials `c1` / `c2`.
|
Chris@87
|
374 The arguments are sequences of coefficients, from lowest order term
|
Chris@87
|
375 to highest, e.g., [1,2,3] represents ``1 + 2*x + 3*x**2``.
|
Chris@87
|
376
|
Chris@87
|
377 Parameters
|
Chris@87
|
378 ----------
|
Chris@87
|
379 c1, c2 : array_like
|
Chris@87
|
380 1-D arrays of polynomial coefficients ordered from low to high.
|
Chris@87
|
381
|
Chris@87
|
382 Returns
|
Chris@87
|
383 -------
|
Chris@87
|
384 [quo, rem] : ndarrays
|
Chris@87
|
385 Of coefficient series representing the quotient and remainder.
|
Chris@87
|
386
|
Chris@87
|
387 See Also
|
Chris@87
|
388 --------
|
Chris@87
|
389 polyadd, polysub, polymul, polypow
|
Chris@87
|
390
|
Chris@87
|
391 Examples
|
Chris@87
|
392 --------
|
Chris@87
|
393 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
394 >>> c1 = (1,2,3)
|
Chris@87
|
395 >>> c2 = (3,2,1)
|
Chris@87
|
396 >>> P.polydiv(c1,c2)
|
Chris@87
|
397 (array([ 3.]), array([-8., -4.]))
|
Chris@87
|
398 >>> P.polydiv(c2,c1)
|
Chris@87
|
399 (array([ 0.33333333]), array([ 2.66666667, 1.33333333]))
|
Chris@87
|
400
|
Chris@87
|
401 """
|
Chris@87
|
402 # c1, c2 are trimmed copies
|
Chris@87
|
403 [c1, c2] = pu.as_series([c1, c2])
|
Chris@87
|
404 if c2[-1] == 0:
|
Chris@87
|
405 raise ZeroDivisionError()
|
Chris@87
|
406
|
Chris@87
|
407 len1 = len(c1)
|
Chris@87
|
408 len2 = len(c2)
|
Chris@87
|
409 if len2 == 1:
|
Chris@87
|
410 return c1/c2[-1], c1[:1]*0
|
Chris@87
|
411 elif len1 < len2:
|
Chris@87
|
412 return c1[:1]*0, c1
|
Chris@87
|
413 else:
|
Chris@87
|
414 dlen = len1 - len2
|
Chris@87
|
415 scl = c2[-1]
|
Chris@87
|
416 c2 = c2[:-1]/scl
|
Chris@87
|
417 i = dlen
|
Chris@87
|
418 j = len1 - 1
|
Chris@87
|
419 while i >= 0:
|
Chris@87
|
420 c1[i:j] -= c2*c1[j]
|
Chris@87
|
421 i -= 1
|
Chris@87
|
422 j -= 1
|
Chris@87
|
423 return c1[j+1:]/scl, pu.trimseq(c1[:j+1])
|
Chris@87
|
424
|
Chris@87
|
425
|
Chris@87
|
426 def polypow(c, pow, maxpower=None):
|
Chris@87
|
427 """Raise a polynomial to a power.
|
Chris@87
|
428
|
Chris@87
|
429 Returns the polynomial `c` raised to the power `pow`. The argument
|
Chris@87
|
430 `c` is a sequence of coefficients ordered from low to high. i.e.,
|
Chris@87
|
431 [1,2,3] is the series ``1 + 2*x + 3*x**2.``
|
Chris@87
|
432
|
Chris@87
|
433 Parameters
|
Chris@87
|
434 ----------
|
Chris@87
|
435 c : array_like
|
Chris@87
|
436 1-D array of array of series coefficients ordered from low to
|
Chris@87
|
437 high degree.
|
Chris@87
|
438 pow : integer
|
Chris@87
|
439 Power to which the series will be raised
|
Chris@87
|
440 maxpower : integer, optional
|
Chris@87
|
441 Maximum power allowed. This is mainly to limit growth of the series
|
Chris@87
|
442 to unmanageable size. Default is 16
|
Chris@87
|
443
|
Chris@87
|
444 Returns
|
Chris@87
|
445 -------
|
Chris@87
|
446 coef : ndarray
|
Chris@87
|
447 Power series of power.
|
Chris@87
|
448
|
Chris@87
|
449 See Also
|
Chris@87
|
450 --------
|
Chris@87
|
451 polyadd, polysub, polymul, polydiv
|
Chris@87
|
452
|
Chris@87
|
453 Examples
|
Chris@87
|
454 --------
|
Chris@87
|
455
|
Chris@87
|
456 """
|
Chris@87
|
457 # c is a trimmed copy
|
Chris@87
|
458 [c] = pu.as_series([c])
|
Chris@87
|
459 power = int(pow)
|
Chris@87
|
460 if power != pow or power < 0:
|
Chris@87
|
461 raise ValueError("Power must be a non-negative integer.")
|
Chris@87
|
462 elif maxpower is not None and power > maxpower:
|
Chris@87
|
463 raise ValueError("Power is too large")
|
Chris@87
|
464 elif power == 0:
|
Chris@87
|
465 return np.array([1], dtype=c.dtype)
|
Chris@87
|
466 elif power == 1:
|
Chris@87
|
467 return c
|
Chris@87
|
468 else:
|
Chris@87
|
469 # This can be made more efficient by using powers of two
|
Chris@87
|
470 # in the usual way.
|
Chris@87
|
471 prd = c
|
Chris@87
|
472 for i in range(2, power + 1):
|
Chris@87
|
473 prd = np.convolve(prd, c)
|
Chris@87
|
474 return prd
|
Chris@87
|
475
|
Chris@87
|
476
|
Chris@87
|
477 def polyder(c, m=1, scl=1, axis=0):
|
Chris@87
|
478 """
|
Chris@87
|
479 Differentiate a polynomial.
|
Chris@87
|
480
|
Chris@87
|
481 Returns the polynomial coefficients `c` differentiated `m` times along
|
Chris@87
|
482 `axis`. At each iteration the result is multiplied by `scl` (the
|
Chris@87
|
483 scaling factor is for use in a linear change of variable). The
|
Chris@87
|
484 argument `c` is an array of coefficients from low to high degree along
|
Chris@87
|
485 each axis, e.g., [1,2,3] represents the polynomial ``1 + 2*x + 3*x**2``
|
Chris@87
|
486 while [[1,2],[1,2]] represents ``1 + 1*x + 2*y + 2*x*y`` if axis=0 is
|
Chris@87
|
487 ``x`` and axis=1 is ``y``.
|
Chris@87
|
488
|
Chris@87
|
489 Parameters
|
Chris@87
|
490 ----------
|
Chris@87
|
491 c : array_like
|
Chris@87
|
492 Array of polynomial coefficients. If c is multidimensional the
|
Chris@87
|
493 different axis correspond to different variables with the degree
|
Chris@87
|
494 in each axis given by the corresponding index.
|
Chris@87
|
495 m : int, optional
|
Chris@87
|
496 Number of derivatives taken, must be non-negative. (Default: 1)
|
Chris@87
|
497 scl : scalar, optional
|
Chris@87
|
498 Each differentiation is multiplied by `scl`. The end result is
|
Chris@87
|
499 multiplication by ``scl**m``. This is for use in a linear change
|
Chris@87
|
500 of variable. (Default: 1)
|
Chris@87
|
501 axis : int, optional
|
Chris@87
|
502 Axis over which the derivative is taken. (Default: 0).
|
Chris@87
|
503
|
Chris@87
|
504 .. versionadded:: 1.7.0
|
Chris@87
|
505
|
Chris@87
|
506 Returns
|
Chris@87
|
507 -------
|
Chris@87
|
508 der : ndarray
|
Chris@87
|
509 Polynomial coefficients of the derivative.
|
Chris@87
|
510
|
Chris@87
|
511 See Also
|
Chris@87
|
512 --------
|
Chris@87
|
513 polyint
|
Chris@87
|
514
|
Chris@87
|
515 Examples
|
Chris@87
|
516 --------
|
Chris@87
|
517 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
518 >>> c = (1,2,3,4) # 1 + 2x + 3x**2 + 4x**3
|
Chris@87
|
519 >>> P.polyder(c) # (d/dx)(c) = 2 + 6x + 12x**2
|
Chris@87
|
520 array([ 2., 6., 12.])
|
Chris@87
|
521 >>> P.polyder(c,3) # (d**3/dx**3)(c) = 24
|
Chris@87
|
522 array([ 24.])
|
Chris@87
|
523 >>> P.polyder(c,scl=-1) # (d/d(-x))(c) = -2 - 6x - 12x**2
|
Chris@87
|
524 array([ -2., -6., -12.])
|
Chris@87
|
525 >>> P.polyder(c,2,-1) # (d**2/d(-x)**2)(c) = 6 + 24x
|
Chris@87
|
526 array([ 6., 24.])
|
Chris@87
|
527
|
Chris@87
|
528 """
|
Chris@87
|
529 c = np.array(c, ndmin=1, copy=1)
|
Chris@87
|
530 if c.dtype.char in '?bBhHiIlLqQpP':
|
Chris@87
|
531 # astype fails with NA
|
Chris@87
|
532 c = c + 0.0
|
Chris@87
|
533 cdt = c.dtype
|
Chris@87
|
534 cnt, iaxis = [int(t) for t in [m, axis]]
|
Chris@87
|
535
|
Chris@87
|
536 if cnt != m:
|
Chris@87
|
537 raise ValueError("The order of derivation must be integer")
|
Chris@87
|
538 if cnt < 0:
|
Chris@87
|
539 raise ValueError("The order of derivation must be non-negative")
|
Chris@87
|
540 if iaxis != axis:
|
Chris@87
|
541 raise ValueError("The axis must be integer")
|
Chris@87
|
542 if not -c.ndim <= iaxis < c.ndim:
|
Chris@87
|
543 raise ValueError("The axis is out of range")
|
Chris@87
|
544 if iaxis < 0:
|
Chris@87
|
545 iaxis += c.ndim
|
Chris@87
|
546
|
Chris@87
|
547 if cnt == 0:
|
Chris@87
|
548 return c
|
Chris@87
|
549
|
Chris@87
|
550 c = np.rollaxis(c, iaxis)
|
Chris@87
|
551 n = len(c)
|
Chris@87
|
552 if cnt >= n:
|
Chris@87
|
553 c = c[:1]*0
|
Chris@87
|
554 else:
|
Chris@87
|
555 for i in range(cnt):
|
Chris@87
|
556 n = n - 1
|
Chris@87
|
557 c *= scl
|
Chris@87
|
558 der = np.empty((n,) + c.shape[1:], dtype=cdt)
|
Chris@87
|
559 for j in range(n, 0, -1):
|
Chris@87
|
560 der[j - 1] = j*c[j]
|
Chris@87
|
561 c = der
|
Chris@87
|
562 c = np.rollaxis(c, 0, iaxis + 1)
|
Chris@87
|
563 return c
|
Chris@87
|
564
|
Chris@87
|
565
|
Chris@87
|
566 def polyint(c, m=1, k=[], lbnd=0, scl=1, axis=0):
|
Chris@87
|
567 """
|
Chris@87
|
568 Integrate a polynomial.
|
Chris@87
|
569
|
Chris@87
|
570 Returns the polynomial coefficients `c` integrated `m` times from
|
Chris@87
|
571 `lbnd` along `axis`. At each iteration the resulting series is
|
Chris@87
|
572 **multiplied** by `scl` and an integration constant, `k`, is added.
|
Chris@87
|
573 The scaling factor is for use in a linear change of variable. ("Buyer
|
Chris@87
|
574 beware": note that, depending on what one is doing, one may want `scl`
|
Chris@87
|
575 to be the reciprocal of what one might expect; for more information,
|
Chris@87
|
576 see the Notes section below.) The argument `c` is an array of
|
Chris@87
|
577 coefficients, from low to high degree along each axis, e.g., [1,2,3]
|
Chris@87
|
578 represents the polynomial ``1 + 2*x + 3*x**2`` while [[1,2],[1,2]]
|
Chris@87
|
579 represents ``1 + 1*x + 2*y + 2*x*y`` if axis=0 is ``x`` and axis=1 is
|
Chris@87
|
580 ``y``.
|
Chris@87
|
581
|
Chris@87
|
582 Parameters
|
Chris@87
|
583 ----------
|
Chris@87
|
584 c : array_like
|
Chris@87
|
585 1-D array of polynomial coefficients, ordered from low to high.
|
Chris@87
|
586 m : int, optional
|
Chris@87
|
587 Order of integration, must be positive. (Default: 1)
|
Chris@87
|
588 k : {[], list, scalar}, optional
|
Chris@87
|
589 Integration constant(s). The value of the first integral at zero
|
Chris@87
|
590 is the first value in the list, the value of the second integral
|
Chris@87
|
591 at zero is the second value, etc. If ``k == []`` (the default),
|
Chris@87
|
592 all constants are set to zero. If ``m == 1``, a single scalar can
|
Chris@87
|
593 be given instead of a list.
|
Chris@87
|
594 lbnd : scalar, optional
|
Chris@87
|
595 The lower bound of the integral. (Default: 0)
|
Chris@87
|
596 scl : scalar, optional
|
Chris@87
|
597 Following each integration the result is *multiplied* by `scl`
|
Chris@87
|
598 before the integration constant is added. (Default: 1)
|
Chris@87
|
599 axis : int, optional
|
Chris@87
|
600 Axis over which the integral is taken. (Default: 0).
|
Chris@87
|
601
|
Chris@87
|
602 .. versionadded:: 1.7.0
|
Chris@87
|
603
|
Chris@87
|
604 Returns
|
Chris@87
|
605 -------
|
Chris@87
|
606 S : ndarray
|
Chris@87
|
607 Coefficient array of the integral.
|
Chris@87
|
608
|
Chris@87
|
609 Raises
|
Chris@87
|
610 ------
|
Chris@87
|
611 ValueError
|
Chris@87
|
612 If ``m < 1``, ``len(k) > m``.
|
Chris@87
|
613
|
Chris@87
|
614 See Also
|
Chris@87
|
615 --------
|
Chris@87
|
616 polyder
|
Chris@87
|
617
|
Chris@87
|
618 Notes
|
Chris@87
|
619 -----
|
Chris@87
|
620 Note that the result of each integration is *multiplied* by `scl`. Why
|
Chris@87
|
621 is this important to note? Say one is making a linear change of
|
Chris@87
|
622 variable :math:`u = ax + b` in an integral relative to `x`. Then
|
Chris@87
|
623 .. math::`dx = du/a`, so one will need to set `scl` equal to
|
Chris@87
|
624 :math:`1/a` - perhaps not what one would have first thought.
|
Chris@87
|
625
|
Chris@87
|
626 Examples
|
Chris@87
|
627 --------
|
Chris@87
|
628 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
629 >>> c = (1,2,3)
|
Chris@87
|
630 >>> P.polyint(c) # should return array([0, 1, 1, 1])
|
Chris@87
|
631 array([ 0., 1., 1., 1.])
|
Chris@87
|
632 >>> P.polyint(c,3) # should return array([0, 0, 0, 1/6, 1/12, 1/20])
|
Chris@87
|
633 array([ 0. , 0. , 0. , 0.16666667, 0.08333333,
|
Chris@87
|
634 0.05 ])
|
Chris@87
|
635 >>> P.polyint(c,k=3) # should return array([3, 1, 1, 1])
|
Chris@87
|
636 array([ 3., 1., 1., 1.])
|
Chris@87
|
637 >>> P.polyint(c,lbnd=-2) # should return array([6, 1, 1, 1])
|
Chris@87
|
638 array([ 6., 1., 1., 1.])
|
Chris@87
|
639 >>> P.polyint(c,scl=-2) # should return array([0, -2, -2, -2])
|
Chris@87
|
640 array([ 0., -2., -2., -2.])
|
Chris@87
|
641
|
Chris@87
|
642 """
|
Chris@87
|
643 c = np.array(c, ndmin=1, copy=1)
|
Chris@87
|
644 if c.dtype.char in '?bBhHiIlLqQpP':
|
Chris@87
|
645 # astype doesn't preserve mask attribute.
|
Chris@87
|
646 c = c + 0.0
|
Chris@87
|
647 cdt = c.dtype
|
Chris@87
|
648 if not np.iterable(k):
|
Chris@87
|
649 k = [k]
|
Chris@87
|
650 cnt, iaxis = [int(t) for t in [m, axis]]
|
Chris@87
|
651
|
Chris@87
|
652 if cnt != m:
|
Chris@87
|
653 raise ValueError("The order of integration must be integer")
|
Chris@87
|
654 if cnt < 0:
|
Chris@87
|
655 raise ValueError("The order of integration must be non-negative")
|
Chris@87
|
656 if len(k) > cnt:
|
Chris@87
|
657 raise ValueError("Too many integration constants")
|
Chris@87
|
658 if iaxis != axis:
|
Chris@87
|
659 raise ValueError("The axis must be integer")
|
Chris@87
|
660 if not -c.ndim <= iaxis < c.ndim:
|
Chris@87
|
661 raise ValueError("The axis is out of range")
|
Chris@87
|
662 if iaxis < 0:
|
Chris@87
|
663 iaxis += c.ndim
|
Chris@87
|
664
|
Chris@87
|
665 if cnt == 0:
|
Chris@87
|
666 return c
|
Chris@87
|
667
|
Chris@87
|
668 k = list(k) + [0]*(cnt - len(k))
|
Chris@87
|
669 c = np.rollaxis(c, iaxis)
|
Chris@87
|
670 for i in range(cnt):
|
Chris@87
|
671 n = len(c)
|
Chris@87
|
672 c *= scl
|
Chris@87
|
673 if n == 1 and np.all(c[0] == 0):
|
Chris@87
|
674 c[0] += k[i]
|
Chris@87
|
675 else:
|
Chris@87
|
676 tmp = np.empty((n + 1,) + c.shape[1:], dtype=cdt)
|
Chris@87
|
677 tmp[0] = c[0]*0
|
Chris@87
|
678 tmp[1] = c[0]
|
Chris@87
|
679 for j in range(1, n):
|
Chris@87
|
680 tmp[j + 1] = c[j]/(j + 1)
|
Chris@87
|
681 tmp[0] += k[i] - polyval(lbnd, tmp)
|
Chris@87
|
682 c = tmp
|
Chris@87
|
683 c = np.rollaxis(c, 0, iaxis + 1)
|
Chris@87
|
684 return c
|
Chris@87
|
685
|
Chris@87
|
686
|
Chris@87
|
687 def polyval(x, c, tensor=True):
|
Chris@87
|
688 """
|
Chris@87
|
689 Evaluate a polynomial at points x.
|
Chris@87
|
690
|
Chris@87
|
691 If `c` is of length `n + 1`, this function returns the value
|
Chris@87
|
692
|
Chris@87
|
693 .. math:: p(x) = c_0 + c_1 * x + ... + c_n * x^n
|
Chris@87
|
694
|
Chris@87
|
695 The parameter `x` is converted to an array only if it is a tuple or a
|
Chris@87
|
696 list, otherwise it is treated as a scalar. In either case, either `x`
|
Chris@87
|
697 or its elements must support multiplication and addition both with
|
Chris@87
|
698 themselves and with the elements of `c`.
|
Chris@87
|
699
|
Chris@87
|
700 If `c` is a 1-D array, then `p(x)` will have the same shape as `x`. If
|
Chris@87
|
701 `c` is multidimensional, then the shape of the result depends on the
|
Chris@87
|
702 value of `tensor`. If `tensor` is true the shape will be c.shape[1:] +
|
Chris@87
|
703 x.shape. If `tensor` is false the shape will be c.shape[1:]. Note that
|
Chris@87
|
704 scalars have shape (,).
|
Chris@87
|
705
|
Chris@87
|
706 Trailing zeros in the coefficients will be used in the evaluation, so
|
Chris@87
|
707 they should be avoided if efficiency is a concern.
|
Chris@87
|
708
|
Chris@87
|
709 Parameters
|
Chris@87
|
710 ----------
|
Chris@87
|
711 x : array_like, compatible object
|
Chris@87
|
712 If `x` is a list or tuple, it is converted to an ndarray, otherwise
|
Chris@87
|
713 it is left unchanged and treated as a scalar. In either case, `x`
|
Chris@87
|
714 or its elements must support addition and multiplication with
|
Chris@87
|
715 with themselves and with the elements of `c`.
|
Chris@87
|
716 c : array_like
|
Chris@87
|
717 Array of coefficients ordered so that the coefficients for terms of
|
Chris@87
|
718 degree n are contained in c[n]. If `c` is multidimensional the
|
Chris@87
|
719 remaining indices enumerate multiple polynomials. In the two
|
Chris@87
|
720 dimensional case the coefficients may be thought of as stored in
|
Chris@87
|
721 the columns of `c`.
|
Chris@87
|
722 tensor : boolean, optional
|
Chris@87
|
723 If True, the shape of the coefficient array is extended with ones
|
Chris@87
|
724 on the right, one for each dimension of `x`. Scalars have dimension 0
|
Chris@87
|
725 for this action. The result is that every column of coefficients in
|
Chris@87
|
726 `c` is evaluated for every element of `x`. If False, `x` is broadcast
|
Chris@87
|
727 over the columns of `c` for the evaluation. This keyword is useful
|
Chris@87
|
728 when `c` is multidimensional. The default value is True.
|
Chris@87
|
729
|
Chris@87
|
730 .. versionadded:: 1.7.0
|
Chris@87
|
731
|
Chris@87
|
732 Returns
|
Chris@87
|
733 -------
|
Chris@87
|
734 values : ndarray, compatible object
|
Chris@87
|
735 The shape of the returned array is described above.
|
Chris@87
|
736
|
Chris@87
|
737 See Also
|
Chris@87
|
738 --------
|
Chris@87
|
739 polyval2d, polygrid2d, polyval3d, polygrid3d
|
Chris@87
|
740
|
Chris@87
|
741 Notes
|
Chris@87
|
742 -----
|
Chris@87
|
743 The evaluation uses Horner's method.
|
Chris@87
|
744
|
Chris@87
|
745 Examples
|
Chris@87
|
746 --------
|
Chris@87
|
747 >>> from numpy.polynomial.polynomial import polyval
|
Chris@87
|
748 >>> polyval(1, [1,2,3])
|
Chris@87
|
749 6.0
|
Chris@87
|
750 >>> a = np.arange(4).reshape(2,2)
|
Chris@87
|
751 >>> a
|
Chris@87
|
752 array([[0, 1],
|
Chris@87
|
753 [2, 3]])
|
Chris@87
|
754 >>> polyval(a, [1,2,3])
|
Chris@87
|
755 array([[ 1., 6.],
|
Chris@87
|
756 [ 17., 34.]])
|
Chris@87
|
757 >>> coef = np.arange(4).reshape(2,2) # multidimensional coefficients
|
Chris@87
|
758 >>> coef
|
Chris@87
|
759 array([[0, 1],
|
Chris@87
|
760 [2, 3]])
|
Chris@87
|
761 >>> polyval([1,2], coef, tensor=True)
|
Chris@87
|
762 array([[ 2., 4.],
|
Chris@87
|
763 [ 4., 7.]])
|
Chris@87
|
764 >>> polyval([1,2], coef, tensor=False)
|
Chris@87
|
765 array([ 2., 7.])
|
Chris@87
|
766
|
Chris@87
|
767 """
|
Chris@87
|
768 c = np.array(c, ndmin=1, copy=0)
|
Chris@87
|
769 if c.dtype.char in '?bBhHiIlLqQpP':
|
Chris@87
|
770 # astype fails with NA
|
Chris@87
|
771 c = c + 0.0
|
Chris@87
|
772 if isinstance(x, (tuple, list)):
|
Chris@87
|
773 x = np.asarray(x)
|
Chris@87
|
774 if isinstance(x, np.ndarray) and tensor:
|
Chris@87
|
775 c = c.reshape(c.shape + (1,)*x.ndim)
|
Chris@87
|
776
|
Chris@87
|
777 c0 = c[-1] + x*0
|
Chris@87
|
778 for i in range(2, len(c) + 1):
|
Chris@87
|
779 c0 = c[-i] + c0*x
|
Chris@87
|
780 return c0
|
Chris@87
|
781
|
Chris@87
|
782
|
Chris@87
|
783 def polyval2d(x, y, c):
|
Chris@87
|
784 """
|
Chris@87
|
785 Evaluate a 2-D polynomial at points (x, y).
|
Chris@87
|
786
|
Chris@87
|
787 This function returns the value
|
Chris@87
|
788
|
Chris@87
|
789 .. math:: p(x,y) = \\sum_{i,j} c_{i,j} * x^i * y^j
|
Chris@87
|
790
|
Chris@87
|
791 The parameters `x` and `y` are converted to arrays only if they are
|
Chris@87
|
792 tuples or a lists, otherwise they are treated as a scalars and they
|
Chris@87
|
793 must have the same shape after conversion. In either case, either `x`
|
Chris@87
|
794 and `y` or their elements must support multiplication and addition both
|
Chris@87
|
795 with themselves and with the elements of `c`.
|
Chris@87
|
796
|
Chris@87
|
797 If `c` has fewer than two dimensions, ones are implicitly appended to
|
Chris@87
|
798 its shape to make it 2-D. The shape of the result will be c.shape[2:] +
|
Chris@87
|
799 x.shape.
|
Chris@87
|
800
|
Chris@87
|
801 Parameters
|
Chris@87
|
802 ----------
|
Chris@87
|
803 x, y : array_like, compatible objects
|
Chris@87
|
804 The two dimensional series is evaluated at the points `(x, y)`,
|
Chris@87
|
805 where `x` and `y` must have the same shape. If `x` or `y` is a list
|
Chris@87
|
806 or tuple, it is first converted to an ndarray, otherwise it is left
|
Chris@87
|
807 unchanged and, if it isn't an ndarray, it is treated as a scalar.
|
Chris@87
|
808 c : array_like
|
Chris@87
|
809 Array of coefficients ordered so that the coefficient of the term
|
Chris@87
|
810 of multi-degree i,j is contained in `c[i,j]`. If `c` has
|
Chris@87
|
811 dimension greater than two the remaining indices enumerate multiple
|
Chris@87
|
812 sets of coefficients.
|
Chris@87
|
813
|
Chris@87
|
814 Returns
|
Chris@87
|
815 -------
|
Chris@87
|
816 values : ndarray, compatible object
|
Chris@87
|
817 The values of the two dimensional polynomial at points formed with
|
Chris@87
|
818 pairs of corresponding values from `x` and `y`.
|
Chris@87
|
819
|
Chris@87
|
820 See Also
|
Chris@87
|
821 --------
|
Chris@87
|
822 polyval, polygrid2d, polyval3d, polygrid3d
|
Chris@87
|
823
|
Chris@87
|
824 Notes
|
Chris@87
|
825 -----
|
Chris@87
|
826
|
Chris@87
|
827 .. versionadded:: 1.7.0
|
Chris@87
|
828
|
Chris@87
|
829 """
|
Chris@87
|
830 try:
|
Chris@87
|
831 x, y = np.array((x, y), copy=0)
|
Chris@87
|
832 except:
|
Chris@87
|
833 raise ValueError('x, y are incompatible')
|
Chris@87
|
834
|
Chris@87
|
835 c = polyval(x, c)
|
Chris@87
|
836 c = polyval(y, c, tensor=False)
|
Chris@87
|
837 return c
|
Chris@87
|
838
|
Chris@87
|
839
|
Chris@87
|
840 def polygrid2d(x, y, c):
|
Chris@87
|
841 """
|
Chris@87
|
842 Evaluate a 2-D polynomial on the Cartesian product of x and y.
|
Chris@87
|
843
|
Chris@87
|
844 This function returns the values:
|
Chris@87
|
845
|
Chris@87
|
846 .. math:: p(a,b) = \\sum_{i,j} c_{i,j} * a^i * b^j
|
Chris@87
|
847
|
Chris@87
|
848 where the points `(a, b)` consist of all pairs formed by taking
|
Chris@87
|
849 `a` from `x` and `b` from `y`. The resulting points form a grid with
|
Chris@87
|
850 `x` in the first dimension and `y` in the second.
|
Chris@87
|
851
|
Chris@87
|
852 The parameters `x` and `y` are converted to arrays only if they are
|
Chris@87
|
853 tuples or a lists, otherwise they are treated as a scalars. In either
|
Chris@87
|
854 case, either `x` and `y` or their elements must support multiplication
|
Chris@87
|
855 and addition both with themselves and with the elements of `c`.
|
Chris@87
|
856
|
Chris@87
|
857 If `c` has fewer than two dimensions, ones are implicitly appended to
|
Chris@87
|
858 its shape to make it 2-D. The shape of the result will be c.shape[2:] +
|
Chris@87
|
859 x.shape + y.shape.
|
Chris@87
|
860
|
Chris@87
|
861 Parameters
|
Chris@87
|
862 ----------
|
Chris@87
|
863 x, y : array_like, compatible objects
|
Chris@87
|
864 The two dimensional series is evaluated at the points in the
|
Chris@87
|
865 Cartesian product of `x` and `y`. If `x` or `y` is a list or
|
Chris@87
|
866 tuple, it is first converted to an ndarray, otherwise it is left
|
Chris@87
|
867 unchanged and, if it isn't an ndarray, it is treated as a scalar.
|
Chris@87
|
868 c : array_like
|
Chris@87
|
869 Array of coefficients ordered so that the coefficients for terms of
|
Chris@87
|
870 degree i,j are contained in ``c[i,j]``. If `c` has dimension
|
Chris@87
|
871 greater than two the remaining indices enumerate multiple sets of
|
Chris@87
|
872 coefficients.
|
Chris@87
|
873
|
Chris@87
|
874 Returns
|
Chris@87
|
875 -------
|
Chris@87
|
876 values : ndarray, compatible object
|
Chris@87
|
877 The values of the two dimensional polynomial at points in the Cartesian
|
Chris@87
|
878 product of `x` and `y`.
|
Chris@87
|
879
|
Chris@87
|
880 See Also
|
Chris@87
|
881 --------
|
Chris@87
|
882 polyval, polyval2d, polyval3d, polygrid3d
|
Chris@87
|
883
|
Chris@87
|
884 Notes
|
Chris@87
|
885 -----
|
Chris@87
|
886
|
Chris@87
|
887 .. versionadded:: 1.7.0
|
Chris@87
|
888
|
Chris@87
|
889 """
|
Chris@87
|
890 c = polyval(x, c)
|
Chris@87
|
891 c = polyval(y, c)
|
Chris@87
|
892 return c
|
Chris@87
|
893
|
Chris@87
|
894
|
Chris@87
|
895 def polyval3d(x, y, z, c):
|
Chris@87
|
896 """
|
Chris@87
|
897 Evaluate a 3-D polynomial at points (x, y, z).
|
Chris@87
|
898
|
Chris@87
|
899 This function returns the values:
|
Chris@87
|
900
|
Chris@87
|
901 .. math:: p(x,y,z) = \\sum_{i,j,k} c_{i,j,k} * x^i * y^j * z^k
|
Chris@87
|
902
|
Chris@87
|
903 The parameters `x`, `y`, and `z` are converted to arrays only if
|
Chris@87
|
904 they are tuples or a lists, otherwise they are treated as a scalars and
|
Chris@87
|
905 they must have the same shape after conversion. In either case, either
|
Chris@87
|
906 `x`, `y`, and `z` or their elements must support multiplication and
|
Chris@87
|
907 addition both with themselves and with the elements of `c`.
|
Chris@87
|
908
|
Chris@87
|
909 If `c` has fewer than 3 dimensions, ones are implicitly appended to its
|
Chris@87
|
910 shape to make it 3-D. The shape of the result will be c.shape[3:] +
|
Chris@87
|
911 x.shape.
|
Chris@87
|
912
|
Chris@87
|
913 Parameters
|
Chris@87
|
914 ----------
|
Chris@87
|
915 x, y, z : array_like, compatible object
|
Chris@87
|
916 The three dimensional series is evaluated at the points
|
Chris@87
|
917 `(x, y, z)`, where `x`, `y`, and `z` must have the same shape. If
|
Chris@87
|
918 any of `x`, `y`, or `z` is a list or tuple, it is first converted
|
Chris@87
|
919 to an ndarray, otherwise it is left unchanged and if it isn't an
|
Chris@87
|
920 ndarray it is treated as a scalar.
|
Chris@87
|
921 c : array_like
|
Chris@87
|
922 Array of coefficients ordered so that the coefficient of the term of
|
Chris@87
|
923 multi-degree i,j,k is contained in ``c[i,j,k]``. If `c` has dimension
|
Chris@87
|
924 greater than 3 the remaining indices enumerate multiple sets of
|
Chris@87
|
925 coefficients.
|
Chris@87
|
926
|
Chris@87
|
927 Returns
|
Chris@87
|
928 -------
|
Chris@87
|
929 values : ndarray, compatible object
|
Chris@87
|
930 The values of the multidimensional polynomial on points formed with
|
Chris@87
|
931 triples of corresponding values from `x`, `y`, and `z`.
|
Chris@87
|
932
|
Chris@87
|
933 See Also
|
Chris@87
|
934 --------
|
Chris@87
|
935 polyval, polyval2d, polygrid2d, polygrid3d
|
Chris@87
|
936
|
Chris@87
|
937 Notes
|
Chris@87
|
938 -----
|
Chris@87
|
939
|
Chris@87
|
940 .. versionadded:: 1.7.0
|
Chris@87
|
941
|
Chris@87
|
942 """
|
Chris@87
|
943 try:
|
Chris@87
|
944 x, y, z = np.array((x, y, z), copy=0)
|
Chris@87
|
945 except:
|
Chris@87
|
946 raise ValueError('x, y, z are incompatible')
|
Chris@87
|
947
|
Chris@87
|
948 c = polyval(x, c)
|
Chris@87
|
949 c = polyval(y, c, tensor=False)
|
Chris@87
|
950 c = polyval(z, c, tensor=False)
|
Chris@87
|
951 return c
|
Chris@87
|
952
|
Chris@87
|
953
|
Chris@87
|
954 def polygrid3d(x, y, z, c):
|
Chris@87
|
955 """
|
Chris@87
|
956 Evaluate a 3-D polynomial on the Cartesian product of x, y and z.
|
Chris@87
|
957
|
Chris@87
|
958 This function returns the values:
|
Chris@87
|
959
|
Chris@87
|
960 .. math:: p(a,b,c) = \\sum_{i,j,k} c_{i,j,k} * a^i * b^j * c^k
|
Chris@87
|
961
|
Chris@87
|
962 where the points `(a, b, c)` consist of all triples formed by taking
|
Chris@87
|
963 `a` from `x`, `b` from `y`, and `c` from `z`. The resulting points form
|
Chris@87
|
964 a grid with `x` in the first dimension, `y` in the second, and `z` in
|
Chris@87
|
965 the third.
|
Chris@87
|
966
|
Chris@87
|
967 The parameters `x`, `y`, and `z` are converted to arrays only if they
|
Chris@87
|
968 are tuples or a lists, otherwise they are treated as a scalars. In
|
Chris@87
|
969 either case, either `x`, `y`, and `z` or their elements must support
|
Chris@87
|
970 multiplication and addition both with themselves and with the elements
|
Chris@87
|
971 of `c`.
|
Chris@87
|
972
|
Chris@87
|
973 If `c` has fewer than three dimensions, ones are implicitly appended to
|
Chris@87
|
974 its shape to make it 3-D. The shape of the result will be c.shape[3:] +
|
Chris@87
|
975 x.shape + y.shape + z.shape.
|
Chris@87
|
976
|
Chris@87
|
977 Parameters
|
Chris@87
|
978 ----------
|
Chris@87
|
979 x, y, z : array_like, compatible objects
|
Chris@87
|
980 The three dimensional series is evaluated at the points in the
|
Chris@87
|
981 Cartesian product of `x`, `y`, and `z`. If `x`,`y`, or `z` is a
|
Chris@87
|
982 list or tuple, it is first converted to an ndarray, otherwise it is
|
Chris@87
|
983 left unchanged and, if it isn't an ndarray, it is treated as a
|
Chris@87
|
984 scalar.
|
Chris@87
|
985 c : array_like
|
Chris@87
|
986 Array of coefficients ordered so that the coefficients for terms of
|
Chris@87
|
987 degree i,j are contained in ``c[i,j]``. If `c` has dimension
|
Chris@87
|
988 greater than two the remaining indices enumerate multiple sets of
|
Chris@87
|
989 coefficients.
|
Chris@87
|
990
|
Chris@87
|
991 Returns
|
Chris@87
|
992 -------
|
Chris@87
|
993 values : ndarray, compatible object
|
Chris@87
|
994 The values of the two dimensional polynomial at points in the Cartesian
|
Chris@87
|
995 product of `x` and `y`.
|
Chris@87
|
996
|
Chris@87
|
997 See Also
|
Chris@87
|
998 --------
|
Chris@87
|
999 polyval, polyval2d, polygrid2d, polyval3d
|
Chris@87
|
1000
|
Chris@87
|
1001 Notes
|
Chris@87
|
1002 -----
|
Chris@87
|
1003
|
Chris@87
|
1004 .. versionadded:: 1.7.0
|
Chris@87
|
1005
|
Chris@87
|
1006 """
|
Chris@87
|
1007 c = polyval(x, c)
|
Chris@87
|
1008 c = polyval(y, c)
|
Chris@87
|
1009 c = polyval(z, c)
|
Chris@87
|
1010 return c
|
Chris@87
|
1011
|
Chris@87
|
1012
|
Chris@87
|
1013 def polyvander(x, deg):
|
Chris@87
|
1014 """Vandermonde matrix of given degree.
|
Chris@87
|
1015
|
Chris@87
|
1016 Returns the Vandermonde matrix of degree `deg` and sample points
|
Chris@87
|
1017 `x`. The Vandermonde matrix is defined by
|
Chris@87
|
1018
|
Chris@87
|
1019 .. math:: V[..., i] = x^i,
|
Chris@87
|
1020
|
Chris@87
|
1021 where `0 <= i <= deg`. The leading indices of `V` index the elements of
|
Chris@87
|
1022 `x` and the last index is the power of `x`.
|
Chris@87
|
1023
|
Chris@87
|
1024 If `c` is a 1-D array of coefficients of length `n + 1` and `V` is the
|
Chris@87
|
1025 matrix ``V = polyvander(x, n)``, then ``np.dot(V, c)`` and
|
Chris@87
|
1026 ``polyval(x, c)`` are the same up to roundoff. This equivalence is
|
Chris@87
|
1027 useful both for least squares fitting and for the evaluation of a large
|
Chris@87
|
1028 number of polynomials of the same degree and sample points.
|
Chris@87
|
1029
|
Chris@87
|
1030 Parameters
|
Chris@87
|
1031 ----------
|
Chris@87
|
1032 x : array_like
|
Chris@87
|
1033 Array of points. The dtype is converted to float64 or complex128
|
Chris@87
|
1034 depending on whether any of the elements are complex. If `x` is
|
Chris@87
|
1035 scalar it is converted to a 1-D array.
|
Chris@87
|
1036 deg : int
|
Chris@87
|
1037 Degree of the resulting matrix.
|
Chris@87
|
1038
|
Chris@87
|
1039 Returns
|
Chris@87
|
1040 -------
|
Chris@87
|
1041 vander : ndarray.
|
Chris@87
|
1042 The Vandermonde matrix. The shape of the returned matrix is
|
Chris@87
|
1043 ``x.shape + (deg + 1,)``, where the last index is the power of `x`.
|
Chris@87
|
1044 The dtype will be the same as the converted `x`.
|
Chris@87
|
1045
|
Chris@87
|
1046 See Also
|
Chris@87
|
1047 --------
|
Chris@87
|
1048 polyvander2d, polyvander3d
|
Chris@87
|
1049
|
Chris@87
|
1050 """
|
Chris@87
|
1051 ideg = int(deg)
|
Chris@87
|
1052 if ideg != deg:
|
Chris@87
|
1053 raise ValueError("deg must be integer")
|
Chris@87
|
1054 if ideg < 0:
|
Chris@87
|
1055 raise ValueError("deg must be non-negative")
|
Chris@87
|
1056
|
Chris@87
|
1057 x = np.array(x, copy=0, ndmin=1) + 0.0
|
Chris@87
|
1058 dims = (ideg + 1,) + x.shape
|
Chris@87
|
1059 dtyp = x.dtype
|
Chris@87
|
1060 v = np.empty(dims, dtype=dtyp)
|
Chris@87
|
1061 v[0] = x*0 + 1
|
Chris@87
|
1062 if ideg > 0:
|
Chris@87
|
1063 v[1] = x
|
Chris@87
|
1064 for i in range(2, ideg + 1):
|
Chris@87
|
1065 v[i] = v[i-1]*x
|
Chris@87
|
1066 return np.rollaxis(v, 0, v.ndim)
|
Chris@87
|
1067
|
Chris@87
|
1068
|
Chris@87
|
1069 def polyvander2d(x, y, deg):
|
Chris@87
|
1070 """Pseudo-Vandermonde matrix of given degrees.
|
Chris@87
|
1071
|
Chris@87
|
1072 Returns the pseudo-Vandermonde matrix of degrees `deg` and sample
|
Chris@87
|
1073 points `(x, y)`. The pseudo-Vandermonde matrix is defined by
|
Chris@87
|
1074
|
Chris@87
|
1075 .. math:: V[..., deg[1]*i + j] = x^i * y^j,
|
Chris@87
|
1076
|
Chris@87
|
1077 where `0 <= i <= deg[0]` and `0 <= j <= deg[1]`. The leading indices of
|
Chris@87
|
1078 `V` index the points `(x, y)` and the last index encodes the powers of
|
Chris@87
|
1079 `x` and `y`.
|
Chris@87
|
1080
|
Chris@87
|
1081 If ``V = polyvander2d(x, y, [xdeg, ydeg])``, then the columns of `V`
|
Chris@87
|
1082 correspond to the elements of a 2-D coefficient array `c` of shape
|
Chris@87
|
1083 (xdeg + 1, ydeg + 1) in the order
|
Chris@87
|
1084
|
Chris@87
|
1085 .. math:: c_{00}, c_{01}, c_{02} ... , c_{10}, c_{11}, c_{12} ...
|
Chris@87
|
1086
|
Chris@87
|
1087 and ``np.dot(V, c.flat)`` and ``polyval2d(x, y, c)`` will be the same
|
Chris@87
|
1088 up to roundoff. This equivalence is useful both for least squares
|
Chris@87
|
1089 fitting and for the evaluation of a large number of 2-D polynomials
|
Chris@87
|
1090 of the same degrees and sample points.
|
Chris@87
|
1091
|
Chris@87
|
1092 Parameters
|
Chris@87
|
1093 ----------
|
Chris@87
|
1094 x, y : array_like
|
Chris@87
|
1095 Arrays of point coordinates, all of the same shape. The dtypes
|
Chris@87
|
1096 will be converted to either float64 or complex128 depending on
|
Chris@87
|
1097 whether any of the elements are complex. Scalars are converted to
|
Chris@87
|
1098 1-D arrays.
|
Chris@87
|
1099 deg : list of ints
|
Chris@87
|
1100 List of maximum degrees of the form [x_deg, y_deg].
|
Chris@87
|
1101
|
Chris@87
|
1102 Returns
|
Chris@87
|
1103 -------
|
Chris@87
|
1104 vander2d : ndarray
|
Chris@87
|
1105 The shape of the returned matrix is ``x.shape + (order,)``, where
|
Chris@87
|
1106 :math:`order = (deg[0]+1)*(deg([1]+1)`. The dtype will be the same
|
Chris@87
|
1107 as the converted `x` and `y`.
|
Chris@87
|
1108
|
Chris@87
|
1109 See Also
|
Chris@87
|
1110 --------
|
Chris@87
|
1111 polyvander, polyvander3d. polyval2d, polyval3d
|
Chris@87
|
1112
|
Chris@87
|
1113 """
|
Chris@87
|
1114 ideg = [int(d) for d in deg]
|
Chris@87
|
1115 is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
|
Chris@87
|
1116 if is_valid != [1, 1]:
|
Chris@87
|
1117 raise ValueError("degrees must be non-negative integers")
|
Chris@87
|
1118 degx, degy = ideg
|
Chris@87
|
1119 x, y = np.array((x, y), copy=0) + 0.0
|
Chris@87
|
1120
|
Chris@87
|
1121 vx = polyvander(x, degx)
|
Chris@87
|
1122 vy = polyvander(y, degy)
|
Chris@87
|
1123 v = vx[..., None]*vy[..., None,:]
|
Chris@87
|
1124 # einsum bug
|
Chris@87
|
1125 #v = np.einsum("...i,...j->...ij", vx, vy)
|
Chris@87
|
1126 return v.reshape(v.shape[:-2] + (-1,))
|
Chris@87
|
1127
|
Chris@87
|
1128
|
Chris@87
|
1129 def polyvander3d(x, y, z, deg):
|
Chris@87
|
1130 """Pseudo-Vandermonde matrix of given degrees.
|
Chris@87
|
1131
|
Chris@87
|
1132 Returns the pseudo-Vandermonde matrix of degrees `deg` and sample
|
Chris@87
|
1133 points `(x, y, z)`. If `l, m, n` are the given degrees in `x, y, z`,
|
Chris@87
|
1134 then The pseudo-Vandermonde matrix is defined by
|
Chris@87
|
1135
|
Chris@87
|
1136 .. math:: V[..., (m+1)(n+1)i + (n+1)j + k] = x^i * y^j * z^k,
|
Chris@87
|
1137
|
Chris@87
|
1138 where `0 <= i <= l`, `0 <= j <= m`, and `0 <= j <= n`. The leading
|
Chris@87
|
1139 indices of `V` index the points `(x, y, z)` and the last index encodes
|
Chris@87
|
1140 the powers of `x`, `y`, and `z`.
|
Chris@87
|
1141
|
Chris@87
|
1142 If ``V = polyvander3d(x, y, z, [xdeg, ydeg, zdeg])``, then the columns
|
Chris@87
|
1143 of `V` correspond to the elements of a 3-D coefficient array `c` of
|
Chris@87
|
1144 shape (xdeg + 1, ydeg + 1, zdeg + 1) in the order
|
Chris@87
|
1145
|
Chris@87
|
1146 .. math:: c_{000}, c_{001}, c_{002},... , c_{010}, c_{011}, c_{012},...
|
Chris@87
|
1147
|
Chris@87
|
1148 and ``np.dot(V, c.flat)`` and ``polyval3d(x, y, z, c)`` will be the
|
Chris@87
|
1149 same up to roundoff. This equivalence is useful both for least squares
|
Chris@87
|
1150 fitting and for the evaluation of a large number of 3-D polynomials
|
Chris@87
|
1151 of the same degrees and sample points.
|
Chris@87
|
1152
|
Chris@87
|
1153 Parameters
|
Chris@87
|
1154 ----------
|
Chris@87
|
1155 x, y, z : array_like
|
Chris@87
|
1156 Arrays of point coordinates, all of the same shape. The dtypes will
|
Chris@87
|
1157 be converted to either float64 or complex128 depending on whether
|
Chris@87
|
1158 any of the elements are complex. Scalars are converted to 1-D
|
Chris@87
|
1159 arrays.
|
Chris@87
|
1160 deg : list of ints
|
Chris@87
|
1161 List of maximum degrees of the form [x_deg, y_deg, z_deg].
|
Chris@87
|
1162
|
Chris@87
|
1163 Returns
|
Chris@87
|
1164 -------
|
Chris@87
|
1165 vander3d : ndarray
|
Chris@87
|
1166 The shape of the returned matrix is ``x.shape + (order,)``, where
|
Chris@87
|
1167 :math:`order = (deg[0]+1)*(deg([1]+1)*(deg[2]+1)`. The dtype will
|
Chris@87
|
1168 be the same as the converted `x`, `y`, and `z`.
|
Chris@87
|
1169
|
Chris@87
|
1170 See Also
|
Chris@87
|
1171 --------
|
Chris@87
|
1172 polyvander, polyvander3d. polyval2d, polyval3d
|
Chris@87
|
1173
|
Chris@87
|
1174 Notes
|
Chris@87
|
1175 -----
|
Chris@87
|
1176
|
Chris@87
|
1177 .. versionadded:: 1.7.0
|
Chris@87
|
1178
|
Chris@87
|
1179 """
|
Chris@87
|
1180 ideg = [int(d) for d in deg]
|
Chris@87
|
1181 is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
|
Chris@87
|
1182 if is_valid != [1, 1, 1]:
|
Chris@87
|
1183 raise ValueError("degrees must be non-negative integers")
|
Chris@87
|
1184 degx, degy, degz = ideg
|
Chris@87
|
1185 x, y, z = np.array((x, y, z), copy=0) + 0.0
|
Chris@87
|
1186
|
Chris@87
|
1187 vx = polyvander(x, degx)
|
Chris@87
|
1188 vy = polyvander(y, degy)
|
Chris@87
|
1189 vz = polyvander(z, degz)
|
Chris@87
|
1190 v = vx[..., None, None]*vy[..., None,:, None]*vz[..., None, None,:]
|
Chris@87
|
1191 # einsum bug
|
Chris@87
|
1192 #v = np.einsum("...i, ...j, ...k->...ijk", vx, vy, vz)
|
Chris@87
|
1193 return v.reshape(v.shape[:-3] + (-1,))
|
Chris@87
|
1194
|
Chris@87
|
1195
|
Chris@87
|
1196 def polyfit(x, y, deg, rcond=None, full=False, w=None):
|
Chris@87
|
1197 """
|
Chris@87
|
1198 Least-squares fit of a polynomial to data.
|
Chris@87
|
1199
|
Chris@87
|
1200 Return the coefficients of a polynomial of degree `deg` that is the
|
Chris@87
|
1201 least squares fit to the data values `y` given at points `x`. If `y` is
|
Chris@87
|
1202 1-D the returned coefficients will also be 1-D. If `y` is 2-D multiple
|
Chris@87
|
1203 fits are done, one for each column of `y`, and the resulting
|
Chris@87
|
1204 coefficients are stored in the corresponding columns of a 2-D return.
|
Chris@87
|
1205 The fitted polynomial(s) are in the form
|
Chris@87
|
1206
|
Chris@87
|
1207 .. math:: p(x) = c_0 + c_1 * x + ... + c_n * x^n,
|
Chris@87
|
1208
|
Chris@87
|
1209 where `n` is `deg`.
|
Chris@87
|
1210
|
Chris@87
|
1211 Parameters
|
Chris@87
|
1212 ----------
|
Chris@87
|
1213 x : array_like, shape (`M`,)
|
Chris@87
|
1214 x-coordinates of the `M` sample (data) points ``(x[i], y[i])``.
|
Chris@87
|
1215 y : array_like, shape (`M`,) or (`M`, `K`)
|
Chris@87
|
1216 y-coordinates of the sample points. Several sets of sample points
|
Chris@87
|
1217 sharing the same x-coordinates can be (independently) fit with one
|
Chris@87
|
1218 call to `polyfit` by passing in for `y` a 2-D array that contains
|
Chris@87
|
1219 one data set per column.
|
Chris@87
|
1220 deg : int
|
Chris@87
|
1221 Degree of the polynomial(s) to be fit.
|
Chris@87
|
1222 rcond : float, optional
|
Chris@87
|
1223 Relative condition number of the fit. Singular values smaller
|
Chris@87
|
1224 than `rcond`, relative to the largest singular value, will be
|
Chris@87
|
1225 ignored. The default value is ``len(x)*eps``, where `eps` is the
|
Chris@87
|
1226 relative precision of the platform's float type, about 2e-16 in
|
Chris@87
|
1227 most cases.
|
Chris@87
|
1228 full : bool, optional
|
Chris@87
|
1229 Switch determining the nature of the return value. When ``False``
|
Chris@87
|
1230 (the default) just the coefficients are returned; when ``True``,
|
Chris@87
|
1231 diagnostic information from the singular value decomposition (used
|
Chris@87
|
1232 to solve the fit's matrix equation) is also returned.
|
Chris@87
|
1233 w : array_like, shape (`M`,), optional
|
Chris@87
|
1234 Weights. If not None, the contribution of each point
|
Chris@87
|
1235 ``(x[i],y[i])`` to the fit is weighted by `w[i]`. Ideally the
|
Chris@87
|
1236 weights are chosen so that the errors of the products ``w[i]*y[i]``
|
Chris@87
|
1237 all have the same variance. The default value is None.
|
Chris@87
|
1238
|
Chris@87
|
1239 .. versionadded:: 1.5.0
|
Chris@87
|
1240
|
Chris@87
|
1241 Returns
|
Chris@87
|
1242 -------
|
Chris@87
|
1243 coef : ndarray, shape (`deg` + 1,) or (`deg` + 1, `K`)
|
Chris@87
|
1244 Polynomial coefficients ordered from low to high. If `y` was 2-D,
|
Chris@87
|
1245 the coefficients in column `k` of `coef` represent the polynomial
|
Chris@87
|
1246 fit to the data in `y`'s `k`-th column.
|
Chris@87
|
1247
|
Chris@87
|
1248 [residuals, rank, singular_values, rcond] : list
|
Chris@87
|
1249 These values are only returned if `full` = True
|
Chris@87
|
1250
|
Chris@87
|
1251 resid -- sum of squared residuals of the least squares fit
|
Chris@87
|
1252 rank -- the numerical rank of the scaled Vandermonde matrix
|
Chris@87
|
1253 sv -- singular values of the scaled Vandermonde matrix
|
Chris@87
|
1254 rcond -- value of `rcond`.
|
Chris@87
|
1255
|
Chris@87
|
1256 For more details, see `linalg.lstsq`.
|
Chris@87
|
1257
|
Chris@87
|
1258 Raises
|
Chris@87
|
1259 ------
|
Chris@87
|
1260 RankWarning
|
Chris@87
|
1261 Raised if the matrix in the least-squares fit is rank deficient.
|
Chris@87
|
1262 The warning is only raised if `full` == False. The warnings can
|
Chris@87
|
1263 be turned off by:
|
Chris@87
|
1264
|
Chris@87
|
1265 >>> import warnings
|
Chris@87
|
1266 >>> warnings.simplefilter('ignore', RankWarning)
|
Chris@87
|
1267
|
Chris@87
|
1268 See Also
|
Chris@87
|
1269 --------
|
Chris@87
|
1270 chebfit, legfit, lagfit, hermfit, hermefit
|
Chris@87
|
1271 polyval : Evaluates a polynomial.
|
Chris@87
|
1272 polyvander : Vandermonde matrix for powers.
|
Chris@87
|
1273 linalg.lstsq : Computes a least-squares fit from the matrix.
|
Chris@87
|
1274 scipy.interpolate.UnivariateSpline : Computes spline fits.
|
Chris@87
|
1275
|
Chris@87
|
1276 Notes
|
Chris@87
|
1277 -----
|
Chris@87
|
1278 The solution is the coefficients of the polynomial `p` that minimizes
|
Chris@87
|
1279 the sum of the weighted squared errors
|
Chris@87
|
1280
|
Chris@87
|
1281 .. math :: E = \\sum_j w_j^2 * |y_j - p(x_j)|^2,
|
Chris@87
|
1282
|
Chris@87
|
1283 where the :math:`w_j` are the weights. This problem is solved by
|
Chris@87
|
1284 setting up the (typically) over-determined matrix equation:
|
Chris@87
|
1285
|
Chris@87
|
1286 .. math :: V(x) * c = w * y,
|
Chris@87
|
1287
|
Chris@87
|
1288 where `V` is the weighted pseudo Vandermonde matrix of `x`, `c` are the
|
Chris@87
|
1289 coefficients to be solved for, `w` are the weights, and `y` are the
|
Chris@87
|
1290 observed values. This equation is then solved using the singular value
|
Chris@87
|
1291 decomposition of `V`.
|
Chris@87
|
1292
|
Chris@87
|
1293 If some of the singular values of `V` are so small that they are
|
Chris@87
|
1294 neglected (and `full` == ``False``), a `RankWarning` will be raised.
|
Chris@87
|
1295 This means that the coefficient values may be poorly determined.
|
Chris@87
|
1296 Fitting to a lower order polynomial will usually get rid of the warning
|
Chris@87
|
1297 (but may not be what you want, of course; if you have independent
|
Chris@87
|
1298 reason(s) for choosing the degree which isn't working, you may have to:
|
Chris@87
|
1299 a) reconsider those reasons, and/or b) reconsider the quality of your
|
Chris@87
|
1300 data). The `rcond` parameter can also be set to a value smaller than
|
Chris@87
|
1301 its default, but the resulting fit may be spurious and have large
|
Chris@87
|
1302 contributions from roundoff error.
|
Chris@87
|
1303
|
Chris@87
|
1304 Polynomial fits using double precision tend to "fail" at about
|
Chris@87
|
1305 (polynomial) degree 20. Fits using Chebyshev or Legendre series are
|
Chris@87
|
1306 generally better conditioned, but much can still depend on the
|
Chris@87
|
1307 distribution of the sample points and the smoothness of the data. If
|
Chris@87
|
1308 the quality of the fit is inadequate, splines may be a good
|
Chris@87
|
1309 alternative.
|
Chris@87
|
1310
|
Chris@87
|
1311 Examples
|
Chris@87
|
1312 --------
|
Chris@87
|
1313 >>> from numpy.polynomial import polynomial as P
|
Chris@87
|
1314 >>> x = np.linspace(-1,1,51) # x "data": [-1, -0.96, ..., 0.96, 1]
|
Chris@87
|
1315 >>> y = x**3 - x + np.random.randn(len(x)) # x^3 - x + N(0,1) "noise"
|
Chris@87
|
1316 >>> c, stats = P.polyfit(x,y,3,full=True)
|
Chris@87
|
1317 >>> c # c[0], c[2] should be approx. 0, c[1] approx. -1, c[3] approx. 1
|
Chris@87
|
1318 array([ 0.01909725, -1.30598256, -0.00577963, 1.02644286])
|
Chris@87
|
1319 >>> stats # note the large SSR, explaining the rather poor results
|
Chris@87
|
1320 [array([ 38.06116253]), 4, array([ 1.38446749, 1.32119158, 0.50443316,
|
Chris@87
|
1321 0.28853036]), 1.1324274851176597e-014]
|
Chris@87
|
1322
|
Chris@87
|
1323 Same thing without the added noise
|
Chris@87
|
1324
|
Chris@87
|
1325 >>> y = x**3 - x
|
Chris@87
|
1326 >>> c, stats = P.polyfit(x,y,3,full=True)
|
Chris@87
|
1327 >>> c # c[0], c[2] should be "very close to 0", c[1] ~= -1, c[3] ~= 1
|
Chris@87
|
1328 array([ -1.73362882e-17, -1.00000000e+00, -2.67471909e-16,
|
Chris@87
|
1329 1.00000000e+00])
|
Chris@87
|
1330 >>> stats # note the minuscule SSR
|
Chris@87
|
1331 [array([ 7.46346754e-31]), 4, array([ 1.38446749, 1.32119158,
|
Chris@87
|
1332 0.50443316, 0.28853036]), 1.1324274851176597e-014]
|
Chris@87
|
1333
|
Chris@87
|
1334 """
|
Chris@87
|
1335 order = int(deg) + 1
|
Chris@87
|
1336 x = np.asarray(x) + 0.0
|
Chris@87
|
1337 y = np.asarray(y) + 0.0
|
Chris@87
|
1338
|
Chris@87
|
1339 # check arguments.
|
Chris@87
|
1340 if deg < 0:
|
Chris@87
|
1341 raise ValueError("expected deg >= 0")
|
Chris@87
|
1342 if x.ndim != 1:
|
Chris@87
|
1343 raise TypeError("expected 1D vector for x")
|
Chris@87
|
1344 if x.size == 0:
|
Chris@87
|
1345 raise TypeError("expected non-empty vector for x")
|
Chris@87
|
1346 if y.ndim < 1 or y.ndim > 2:
|
Chris@87
|
1347 raise TypeError("expected 1D or 2D array for y")
|
Chris@87
|
1348 if len(x) != len(y):
|
Chris@87
|
1349 raise TypeError("expected x and y to have same length")
|
Chris@87
|
1350
|
Chris@87
|
1351 # set up the least squares matrices in transposed form
|
Chris@87
|
1352 lhs = polyvander(x, deg).T
|
Chris@87
|
1353 rhs = y.T
|
Chris@87
|
1354 if w is not None:
|
Chris@87
|
1355 w = np.asarray(w) + 0.0
|
Chris@87
|
1356 if w.ndim != 1:
|
Chris@87
|
1357 raise TypeError("expected 1D vector for w")
|
Chris@87
|
1358 if len(x) != len(w):
|
Chris@87
|
1359 raise TypeError("expected x and w to have same length")
|
Chris@87
|
1360 # apply weights. Don't use inplace operations as they
|
Chris@87
|
1361 # can cause problems with NA.
|
Chris@87
|
1362 lhs = lhs * w
|
Chris@87
|
1363 rhs = rhs * w
|
Chris@87
|
1364
|
Chris@87
|
1365 # set rcond
|
Chris@87
|
1366 if rcond is None:
|
Chris@87
|
1367 rcond = len(x)*np.finfo(x.dtype).eps
|
Chris@87
|
1368
|
Chris@87
|
1369 # Determine the norms of the design matrix columns.
|
Chris@87
|
1370 if issubclass(lhs.dtype.type, np.complexfloating):
|
Chris@87
|
1371 scl = np.sqrt((np.square(lhs.real) + np.square(lhs.imag)).sum(1))
|
Chris@87
|
1372 else:
|
Chris@87
|
1373 scl = np.sqrt(np.square(lhs).sum(1))
|
Chris@87
|
1374 scl[scl == 0] = 1
|
Chris@87
|
1375
|
Chris@87
|
1376 # Solve the least squares problem.
|
Chris@87
|
1377 c, resids, rank, s = la.lstsq(lhs.T/scl, rhs.T, rcond)
|
Chris@87
|
1378 c = (c.T/scl).T
|
Chris@87
|
1379
|
Chris@87
|
1380 # warn on rank reduction
|
Chris@87
|
1381 if rank != order and not full:
|
Chris@87
|
1382 msg = "The fit may be poorly conditioned"
|
Chris@87
|
1383 warnings.warn(msg, pu.RankWarning)
|
Chris@87
|
1384
|
Chris@87
|
1385 if full:
|
Chris@87
|
1386 return c, [resids, rank, s, rcond]
|
Chris@87
|
1387 else:
|
Chris@87
|
1388 return c
|
Chris@87
|
1389
|
Chris@87
|
1390
|
Chris@87
|
1391 def polycompanion(c):
|
Chris@87
|
1392 """
|
Chris@87
|
1393 Return the companion matrix of c.
|
Chris@87
|
1394
|
Chris@87
|
1395 The companion matrix for power series cannot be made symmetric by
|
Chris@87
|
1396 scaling the basis, so this function differs from those for the
|
Chris@87
|
1397 orthogonal polynomials.
|
Chris@87
|
1398
|
Chris@87
|
1399 Parameters
|
Chris@87
|
1400 ----------
|
Chris@87
|
1401 c : array_like
|
Chris@87
|
1402 1-D array of polynomial coefficients ordered from low to high
|
Chris@87
|
1403 degree.
|
Chris@87
|
1404
|
Chris@87
|
1405 Returns
|
Chris@87
|
1406 -------
|
Chris@87
|
1407 mat : ndarray
|
Chris@87
|
1408 Companion matrix of dimensions (deg, deg).
|
Chris@87
|
1409
|
Chris@87
|
1410 Notes
|
Chris@87
|
1411 -----
|
Chris@87
|
1412
|
Chris@87
|
1413 .. versionadded:: 1.7.0
|
Chris@87
|
1414
|
Chris@87
|
1415 """
|
Chris@87
|
1416 # c is a trimmed copy
|
Chris@87
|
1417 [c] = pu.as_series([c])
|
Chris@87
|
1418 if len(c) < 2:
|
Chris@87
|
1419 raise ValueError('Series must have maximum degree of at least 1.')
|
Chris@87
|
1420 if len(c) == 2:
|
Chris@87
|
1421 return np.array([[-c[0]/c[1]]])
|
Chris@87
|
1422
|
Chris@87
|
1423 n = len(c) - 1
|
Chris@87
|
1424 mat = np.zeros((n, n), dtype=c.dtype)
|
Chris@87
|
1425 bot = mat.reshape(-1)[n::n+1]
|
Chris@87
|
1426 bot[...] = 1
|
Chris@87
|
1427 mat[:, -1] -= c[:-1]/c[-1]
|
Chris@87
|
1428 return mat
|
Chris@87
|
1429
|
Chris@87
|
1430
|
Chris@87
|
1431 def polyroots(c):
|
Chris@87
|
1432 """
|
Chris@87
|
1433 Compute the roots of a polynomial.
|
Chris@87
|
1434
|
Chris@87
|
1435 Return the roots (a.k.a. "zeros") of the polynomial
|
Chris@87
|
1436
|
Chris@87
|
1437 .. math:: p(x) = \\sum_i c[i] * x^i.
|
Chris@87
|
1438
|
Chris@87
|
1439 Parameters
|
Chris@87
|
1440 ----------
|
Chris@87
|
1441 c : 1-D array_like
|
Chris@87
|
1442 1-D array of polynomial coefficients.
|
Chris@87
|
1443
|
Chris@87
|
1444 Returns
|
Chris@87
|
1445 -------
|
Chris@87
|
1446 out : ndarray
|
Chris@87
|
1447 Array of the roots of the polynomial. If all the roots are real,
|
Chris@87
|
1448 then `out` is also real, otherwise it is complex.
|
Chris@87
|
1449
|
Chris@87
|
1450 See Also
|
Chris@87
|
1451 --------
|
Chris@87
|
1452 chebroots
|
Chris@87
|
1453
|
Chris@87
|
1454 Notes
|
Chris@87
|
1455 -----
|
Chris@87
|
1456 The root estimates are obtained as the eigenvalues of the companion
|
Chris@87
|
1457 matrix, Roots far from the origin of the complex plane may have large
|
Chris@87
|
1458 errors due to the numerical instability of the power series for such
|
Chris@87
|
1459 values. Roots with multiplicity greater than 1 will also show larger
|
Chris@87
|
1460 errors as the value of the series near such points is relatively
|
Chris@87
|
1461 insensitive to errors in the roots. Isolated roots near the origin can
|
Chris@87
|
1462 be improved by a few iterations of Newton's method.
|
Chris@87
|
1463
|
Chris@87
|
1464 Examples
|
Chris@87
|
1465 --------
|
Chris@87
|
1466 >>> import numpy.polynomial.polynomial as poly
|
Chris@87
|
1467 >>> poly.polyroots(poly.polyfromroots((-1,0,1)))
|
Chris@87
|
1468 array([-1., 0., 1.])
|
Chris@87
|
1469 >>> poly.polyroots(poly.polyfromroots((-1,0,1))).dtype
|
Chris@87
|
1470 dtype('float64')
|
Chris@87
|
1471 >>> j = complex(0,1)
|
Chris@87
|
1472 >>> poly.polyroots(poly.polyfromroots((-j,0,j)))
|
Chris@87
|
1473 array([ 0.00000000e+00+0.j, 0.00000000e+00+1.j, 2.77555756e-17-1.j])
|
Chris@87
|
1474
|
Chris@87
|
1475 """
|
Chris@87
|
1476 # c is a trimmed copy
|
Chris@87
|
1477 [c] = pu.as_series([c])
|
Chris@87
|
1478 if len(c) < 2:
|
Chris@87
|
1479 return np.array([], dtype=c.dtype)
|
Chris@87
|
1480 if len(c) == 2:
|
Chris@87
|
1481 return np.array([-c[0]/c[1]])
|
Chris@87
|
1482
|
Chris@87
|
1483 m = polycompanion(c)
|
Chris@87
|
1484 r = la.eigvals(m)
|
Chris@87
|
1485 r.sort()
|
Chris@87
|
1486 return r
|
Chris@87
|
1487
|
Chris@87
|
1488
|
Chris@87
|
1489 #
|
Chris@87
|
1490 # polynomial class
|
Chris@87
|
1491 #
|
Chris@87
|
1492
|
Chris@87
|
1493 class Polynomial(ABCPolyBase):
|
Chris@87
|
1494 """A power series class.
|
Chris@87
|
1495
|
Chris@87
|
1496 The Polynomial class provides the standard Python numerical methods
|
Chris@87
|
1497 '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
|
Chris@87
|
1498 attributes and methods listed in the `ABCPolyBase` documentation.
|
Chris@87
|
1499
|
Chris@87
|
1500 Parameters
|
Chris@87
|
1501 ----------
|
Chris@87
|
1502 coef : array_like
|
Chris@87
|
1503 Polynomial coefficients in order of increasing degree, i.e.,
|
Chris@87
|
1504 ``(1, 2, 3)`` give ``1 + 2*x + 3*x**2``.
|
Chris@87
|
1505 domain : (2,) array_like, optional
|
Chris@87
|
1506 Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
|
Chris@87
|
1507 to the interval ``[window[0], window[1]]`` by shifting and scaling.
|
Chris@87
|
1508 The default value is [-1, 1].
|
Chris@87
|
1509 window : (2,) array_like, optional
|
Chris@87
|
1510 Window, see `domain` for its use. The default value is [-1, 1].
|
Chris@87
|
1511
|
Chris@87
|
1512 .. versionadded:: 1.6.0
|
Chris@87
|
1513
|
Chris@87
|
1514 """
|
Chris@87
|
1515 # Virtual Functions
|
Chris@87
|
1516 _add = staticmethod(polyadd)
|
Chris@87
|
1517 _sub = staticmethod(polysub)
|
Chris@87
|
1518 _mul = staticmethod(polymul)
|
Chris@87
|
1519 _div = staticmethod(polydiv)
|
Chris@87
|
1520 _pow = staticmethod(polypow)
|
Chris@87
|
1521 _val = staticmethod(polyval)
|
Chris@87
|
1522 _int = staticmethod(polyint)
|
Chris@87
|
1523 _der = staticmethod(polyder)
|
Chris@87
|
1524 _fit = staticmethod(polyfit)
|
Chris@87
|
1525 _line = staticmethod(polyline)
|
Chris@87
|
1526 _roots = staticmethod(polyroots)
|
Chris@87
|
1527 _fromroots = staticmethod(polyfromroots)
|
Chris@87
|
1528
|
Chris@87
|
1529 # Virtual properties
|
Chris@87
|
1530 nickname = 'poly'
|
Chris@87
|
1531 domain = np.array(polydomain)
|
Chris@87
|
1532 window = np.array(polydomain)
|