hermite_e.py 51 KB

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