legendre.py 49 KB

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