polynomial.py 51 KB

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