laguerre.py 51 KB

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