hermite.py 51 KB

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