hermite.py 53 KB

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