hermite_e.py 51 KB

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