hermite_e.py 57 KB

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