legendre.py 56 KB

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