polynomial.py 52 KB

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