legendre.py 50 KB

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