hermite_e.py 51 KB

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