legendre.py 50 KB

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