orthogonal.py 58 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613161416151616161716181619162016211622162316241625162616271628162916301631163216331634163516361637163816391640164116421643164416451646164716481649165016511652165316541655165616571658165916601661166216631664166516661667166816691670167116721673167416751676167716781679168016811682168316841685168616871688168916901691169216931694169516961697169816991700170117021703170417051706170717081709171017111712171317141715171617171718171917201721172217231724172517261727172817291730173117321733173417351736173717381739174017411742174317441745174617471748174917501751175217531754175517561757175817591760176117621763176417651766176717681769177017711772177317741775177617771778177917801781178217831784178517861787178817891790179117921793179417951796179717981799180018011802180318041805180618071808180918101811181218131814181518161817181818191820182118221823182418251826182718281829183018311832183318341835183618371838183918401841184218431844184518461847184818491850185118521853185418551856185718581859186018611862186318641865186618671868186918701871187218731874187518761877187818791880188118821883188418851886188718881889189018911892189318941895189618971898189919001901190219031904190519061907190819091910191119121913191419151916191719181919192019211922192319241925192619271928192919301931193219331934193519361937193819391940194119421943194419451946194719481949195019511952195319541955195619571958195919601961196219631964196519661967196819691970197119721973197419751976197719781979198019811982198319841985198619871988198919901991199219931994199519961997199819992000200120022003200420052006200720082009201020112012201320142015201620172018201920202021202220232024202520262027202820292030203120322033203420352036203720382039204020412042204320442045204620472048204920502051205220532054205520562057205820592060206120622063206420652066206720682069207020712072207320742075207620772078207920802081
  1. """
  2. A collection of functions to find the weights and abscissas for
  3. Gaussian Quadrature.
  4. These calculations are done by finding the eigenvalues of a
  5. tridiagonal matrix whose entries are dependent on the coefficients
  6. in the recursion formula for the orthogonal polynomials with the
  7. corresponding weighting function over the interval.
  8. Many recursion relations for orthogonal polynomials are given:
  9. .. math::
  10. a1n f_{n+1} (x) = (a2n + a3n x ) f_n (x) - a4n f_{n-1} (x)
  11. The recursion relation of interest is
  12. .. math::
  13. P_{n+1} (x) = (x - A_n) P_n (x) - B_n P_{n-1} (x)
  14. where :math:`P` has a different normalization than :math:`f`.
  15. The coefficients can be found as:
  16. .. math::
  17. A_n = -a2n / a3n
  18. \\qquad
  19. B_n = ( a4n / a3n \\sqrt{h_n-1 / h_n})^2
  20. where
  21. .. math::
  22. h_n = \\int_a^b w(x) f_n(x)^2
  23. assume:
  24. .. math::
  25. P_0 (x) = 1
  26. \\qquad
  27. P_{-1} (x) == 0
  28. For the mathematical background, see [golub.welsch-1969-mathcomp]_ and
  29. [abramowitz.stegun-1965]_.
  30. References
  31. ----------
  32. .. [golub.welsch-1969-mathcomp]
  33. Golub, Gene H, and John H Welsch. 1969. Calculation of Gauss
  34. Quadrature Rules. *Mathematics of Computation* 23, 221-230+s1--s10.
  35. .. [abramowitz.stegun-1965]
  36. Abramowitz, Milton, and Irene A Stegun. (1965) *Handbook of
  37. Mathematical Functions: with Formulas, Graphs, and Mathematical
  38. Tables*. Gaithersburg, MD: National Bureau of Standards.
  39. http://www.math.sfu.ca/~cbm/aands/
  40. .. [townsend.trogdon.olver-2014]
  41. Townsend, A. and Trogdon, T. and Olver, S. (2014)
  42. *Fast computation of Gauss quadrature nodes and
  43. weights on the whole real line*. :arXiv:`1410.5286`.
  44. .. [townsend.trogdon.olver-2015]
  45. Townsend, A. and Trogdon, T. and Olver, S. (2015)
  46. *Fast computation of Gauss quadrature nodes and
  47. weights on the whole real line*.
  48. IMA Journal of Numerical Analysis
  49. :doi:`10.1093/imanum/drv002`.
  50. """
  51. #
  52. # Author: Travis Oliphant 2000
  53. # Updated Sep. 2003 (fixed bugs --- tested to be accurate)
  54. from __future__ import division, print_function, absolute_import
  55. # Scipy imports.
  56. import numpy as np
  57. from numpy import (exp, inf, pi, sqrt, floor, sin, cos, around, int,
  58. hstack, arccos, arange)
  59. from scipy import linalg
  60. from scipy.special import airy
  61. # Local imports.
  62. from . import _ufuncs
  63. from . import _ufuncs as cephes
  64. _gam = cephes.gamma
  65. from . import specfun
  66. _polyfuns = ['legendre', 'chebyt', 'chebyu', 'chebyc', 'chebys',
  67. 'jacobi', 'laguerre', 'genlaguerre', 'hermite',
  68. 'hermitenorm', 'gegenbauer', 'sh_legendre', 'sh_chebyt',
  69. 'sh_chebyu', 'sh_jacobi']
  70. # Correspondence between new and old names of root functions
  71. _rootfuns_map = {'roots_legendre': 'p_roots',
  72. 'roots_chebyt': 't_roots',
  73. 'roots_chebyu': 'u_roots',
  74. 'roots_chebyc': 'c_roots',
  75. 'roots_chebys': 's_roots',
  76. 'roots_jacobi': 'j_roots',
  77. 'roots_laguerre': 'l_roots',
  78. 'roots_genlaguerre': 'la_roots',
  79. 'roots_hermite': 'h_roots',
  80. 'roots_hermitenorm': 'he_roots',
  81. 'roots_gegenbauer': 'cg_roots',
  82. 'roots_sh_legendre': 'ps_roots',
  83. 'roots_sh_chebyt': 'ts_roots',
  84. 'roots_sh_chebyu': 'us_roots',
  85. 'roots_sh_jacobi': 'js_roots'}
  86. _evalfuns = ['eval_legendre', 'eval_chebyt', 'eval_chebyu',
  87. 'eval_chebyc', 'eval_chebys', 'eval_jacobi',
  88. 'eval_laguerre', 'eval_genlaguerre', 'eval_hermite',
  89. 'eval_hermitenorm', 'eval_gegenbauer',
  90. 'eval_sh_legendre', 'eval_sh_chebyt', 'eval_sh_chebyu',
  91. 'eval_sh_jacobi']
  92. __all__ = _polyfuns + list(_rootfuns_map.keys()) + _evalfuns + ['poch', 'binom']
  93. class orthopoly1d(np.poly1d):
  94. def __init__(self, roots, weights=None, hn=1.0, kn=1.0, wfunc=None,
  95. limits=None, monic=False, eval_func=None):
  96. equiv_weights = [weights[k] / wfunc(roots[k]) for
  97. k in range(len(roots))]
  98. mu = sqrt(hn)
  99. if monic:
  100. evf = eval_func
  101. if evf:
  102. knn = kn
  103. eval_func = lambda x: evf(x) / knn
  104. mu = mu / abs(kn)
  105. kn = 1.0
  106. # compute coefficients from roots, then scale
  107. poly = np.poly1d(roots, r=True)
  108. np.poly1d.__init__(self, poly.coeffs * float(kn))
  109. # TODO: In numpy 1.13, there is no need to use __dict__ to access attributes
  110. self.__dict__['weights'] = np.array(list(zip(roots,
  111. weights, equiv_weights)))
  112. self.__dict__['weight_func'] = wfunc
  113. self.__dict__['limits'] = limits
  114. self.__dict__['normcoef'] = mu
  115. # Note: eval_func will be discarded on arithmetic
  116. self.__dict__['_eval_func'] = eval_func
  117. def __call__(self, v):
  118. if self._eval_func and not isinstance(v, np.poly1d):
  119. return self._eval_func(v)
  120. else:
  121. return np.poly1d.__call__(self, v)
  122. def _scale(self, p):
  123. if p == 1.0:
  124. return
  125. try:
  126. self._coeffs
  127. except AttributeError:
  128. self.__dict__['coeffs'] *= p
  129. else:
  130. # the coeffs attr is be made private in future versions of numpy
  131. self._coeffs *= p
  132. evf = self._eval_func
  133. if evf:
  134. self.__dict__['_eval_func'] = lambda x: evf(x) * p
  135. self.__dict__['normcoef'] *= p
  136. def _gen_roots_and_weights(n, mu0, an_func, bn_func, f, df, symmetrize, mu):
  137. """[x,w] = gen_roots_and_weights(n,an_func,sqrt_bn_func,mu)
  138. Returns the roots (x) of an nth order orthogonal polynomial,
  139. and weights (w) to use in appropriate Gaussian quadrature with that
  140. orthogonal polynomial.
  141. The polynomials have the recurrence relation
  142. P_n+1(x) = (x - A_n) P_n(x) - B_n P_n-1(x)
  143. an_func(n) should return A_n
  144. sqrt_bn_func(n) should return sqrt(B_n)
  145. mu ( = h_0 ) is the integral of the weight over the orthogonal
  146. interval
  147. """
  148. k = np.arange(n, dtype='d')
  149. c = np.zeros((2, n))
  150. c[0,1:] = bn_func(k[1:])
  151. c[1,:] = an_func(k)
  152. x = linalg.eigvals_banded(c, overwrite_a_band=True)
  153. # improve roots by one application of Newton's method
  154. y = f(n, x)
  155. dy = df(n, x)
  156. x -= y/dy
  157. fm = f(n-1, x)
  158. fm /= np.abs(fm).max()
  159. dy /= np.abs(dy).max()
  160. w = 1.0 / (fm * dy)
  161. if symmetrize:
  162. w = (w + w[::-1]) / 2
  163. x = (x - x[::-1]) / 2
  164. w *= mu0 / w.sum()
  165. if mu:
  166. return x, w, mu0
  167. else:
  168. return x, w
  169. # Jacobi Polynomials 1 P^(alpha,beta)_n(x)
  170. def roots_jacobi(n, alpha, beta, mu=False):
  171. r"""Gauss-Jacobi quadrature.
  172. Computes the sample points and weights for Gauss-Jacobi quadrature. The
  173. sample points are the roots of the n-th degree Jacobi polynomial,
  174. :math:`P^{\alpha, \beta}_n(x)`. These sample points and weights
  175. correctly integrate polynomials of degree :math:`2n - 1` or less over the
  176. interval :math:`[-1, 1]` with weight function
  177. :math:`f(x) = (1 - x)^{\alpha} (1 + x)^{\beta}`.
  178. Parameters
  179. ----------
  180. n : int
  181. quadrature order
  182. alpha : float
  183. alpha must be > -1
  184. beta : float
  185. beta must be > -1
  186. mu : bool, optional
  187. If True, return the sum of the weights, optional.
  188. Returns
  189. -------
  190. x : ndarray
  191. Sample points
  192. w : ndarray
  193. Weights
  194. mu : float
  195. Sum of the weights
  196. See Also
  197. --------
  198. scipy.integrate.quadrature
  199. scipy.integrate.fixed_quad
  200. """
  201. m = int(n)
  202. if n < 1 or n != m:
  203. raise ValueError("n must be a positive integer.")
  204. if alpha <= -1 or beta <= -1:
  205. raise ValueError("alpha and beta must be greater than -1.")
  206. if alpha == 0.0 and beta == 0.0:
  207. return roots_legendre(m, mu)
  208. if alpha == beta:
  209. return roots_gegenbauer(m, alpha+0.5, mu)
  210. mu0 = 2.0**(alpha+beta+1)*cephes.beta(alpha+1, beta+1)
  211. a = alpha
  212. b = beta
  213. if a + b == 0.0:
  214. an_func = lambda k: np.where(k == 0, (b-a)/(2+a+b), 0.0)
  215. else:
  216. an_func = lambda k: np.where(k == 0, (b-a)/(2+a+b),
  217. (b*b - a*a) / ((2.0*k+a+b)*(2.0*k+a+b+2)))
  218. bn_func = lambda k: 2.0 / (2.0*k+a+b)*np.sqrt((k+a)*(k+b) / (2*k+a+b+1)) \
  219. * np.where(k == 1, 1.0, np.sqrt(k*(k+a+b) / (2.0*k+a+b-1)))
  220. f = lambda n, x: cephes.eval_jacobi(n, a, b, x)
  221. df = lambda n, x: 0.5 * (n + a + b + 1) \
  222. * cephes.eval_jacobi(n-1, a+1, b+1, x)
  223. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, False, mu)
  224. def jacobi(n, alpha, beta, monic=False):
  225. r"""Jacobi polynomial.
  226. Defined to be the solution of
  227. .. math::
  228. (1 - x^2)\frac{d^2}{dx^2}P_n^{(\alpha, \beta)}
  229. + (\beta - \alpha - (\alpha + \beta + 2)x)
  230. \frac{d}{dx}P_n^{(\alpha, \beta)}
  231. + n(n + \alpha + \beta + 1)P_n^{(\alpha, \beta)} = 0
  232. for :math:`\alpha, \beta > -1`; :math:`P_n^{(\alpha, \beta)}` is a
  233. polynomial of degree :math:`n`.
  234. Parameters
  235. ----------
  236. n : int
  237. Degree of the polynomial.
  238. alpha : float
  239. Parameter, must be greater than -1.
  240. beta : float
  241. Parameter, must be greater than -1.
  242. monic : bool, optional
  243. If `True`, scale the leading coefficient to be 1. Default is
  244. `False`.
  245. Returns
  246. -------
  247. P : orthopoly1d
  248. Jacobi polynomial.
  249. Notes
  250. -----
  251. For fixed :math:`\alpha, \beta`, the polynomials
  252. :math:`P_n^{(\alpha, \beta)}` are orthogonal over :math:`[-1, 1]`
  253. with weight function :math:`(1 - x)^\alpha(1 + x)^\beta`.
  254. """
  255. if n < 0:
  256. raise ValueError("n must be nonnegative.")
  257. wfunc = lambda x: (1 - x)**alpha * (1 + x)**beta
  258. if n == 0:
  259. return orthopoly1d([], [], 1.0, 1.0, wfunc, (-1, 1), monic,
  260. eval_func=np.ones_like)
  261. x, w, mu = roots_jacobi(n, alpha, beta, mu=True)
  262. ab1 = alpha + beta + 1.0
  263. hn = 2**ab1 / (2 * n + ab1) * _gam(n + alpha + 1)
  264. hn *= _gam(n + beta + 1.0) / _gam(n + 1) / _gam(n + ab1)
  265. kn = _gam(2 * n + ab1) / 2.0**n / _gam(n + 1) / _gam(n + ab1)
  266. # here kn = coefficient on x^n term
  267. p = orthopoly1d(x, w, hn, kn, wfunc, (-1, 1), monic,
  268. lambda x: eval_jacobi(n, alpha, beta, x))
  269. return p
  270. # Jacobi Polynomials shifted G_n(p,q,x)
  271. def roots_sh_jacobi(n, p1, q1, mu=False):
  272. """Gauss-Jacobi (shifted) quadrature.
  273. Computes the sample points and weights for Gauss-Jacobi (shifted)
  274. quadrature. The sample points are the roots of the n-th degree shifted
  275. Jacobi polynomial, :math:`G^{p,q}_n(x)`. These sample points and weights
  276. correctly integrate polynomials of degree :math:`2n - 1` or less over the
  277. interval :math:`[0, 1]` with weight function
  278. :math:`f(x) = (1 - x)^{p-q} x^{q-1}`
  279. Parameters
  280. ----------
  281. n : int
  282. quadrature order
  283. p1 : float
  284. (p1 - q1) must be > -1
  285. q1 : float
  286. q1 must be > 0
  287. mu : bool, optional
  288. If True, return the sum of the weights, optional.
  289. Returns
  290. -------
  291. x : ndarray
  292. Sample points
  293. w : ndarray
  294. Weights
  295. mu : float
  296. Sum of the weights
  297. See Also
  298. --------
  299. scipy.integrate.quadrature
  300. scipy.integrate.fixed_quad
  301. """
  302. if (p1-q1) <= -1 or q1 <= 0:
  303. raise ValueError("(p - q) must be greater than -1, and q must be greater than 0.")
  304. x, w, m = roots_jacobi(n, p1-q1, q1-1, True)
  305. x = (x + 1) / 2
  306. scale = 2.0**p1
  307. w /= scale
  308. m /= scale
  309. if mu:
  310. return x, w, m
  311. else:
  312. return x, w
  313. def sh_jacobi(n, p, q, monic=False):
  314. r"""Shifted Jacobi polynomial.
  315. Defined by
  316. .. math::
  317. G_n^{(p, q)}(x)
  318. = \binom{2n + p - 1}{n}^{-1}P_n^{(p - q, q - 1)}(2x - 1),
  319. where :math:`P_n^{(\cdot, \cdot)}` is the nth Jacobi polynomial.
  320. Parameters
  321. ----------
  322. n : int
  323. Degree of the polynomial.
  324. p : float
  325. Parameter, must have :math:`p > q - 1`.
  326. q : float
  327. Parameter, must be greater than 0.
  328. monic : bool, optional
  329. If `True`, scale the leading coefficient to be 1. Default is
  330. `False`.
  331. Returns
  332. -------
  333. G : orthopoly1d
  334. Shifted Jacobi polynomial.
  335. Notes
  336. -----
  337. For fixed :math:`p, q`, the polynomials :math:`G_n^{(p, q)}` are
  338. orthogonal over :math:`[0, 1]` with weight function :math:`(1 -
  339. x)^{p - q}x^{q - 1}`.
  340. """
  341. if n < 0:
  342. raise ValueError("n must be nonnegative.")
  343. wfunc = lambda x: (1.0 - x)**(p - q) * (x)**(q - 1.)
  344. if n == 0:
  345. return orthopoly1d([], [], 1.0, 1.0, wfunc, (-1, 1), monic,
  346. eval_func=np.ones_like)
  347. n1 = n
  348. x, w, mu0 = roots_sh_jacobi(n1, p, q, mu=True)
  349. hn = _gam(n + 1) * _gam(n + q) * _gam(n + p) * _gam(n + p - q + 1)
  350. hn /= (2 * n + p) * (_gam(2 * n + p)**2)
  351. # kn = 1.0 in standard form so monic is redundant. Kept for compatibility.
  352. kn = 1.0
  353. pp = orthopoly1d(x, w, hn, kn, wfunc=wfunc, limits=(0, 1), monic=monic,
  354. eval_func=lambda x: eval_sh_jacobi(n, p, q, x))
  355. return pp
  356. # Generalized Laguerre L^(alpha)_n(x)
  357. def roots_genlaguerre(n, alpha, mu=False):
  358. r"""Gauss-generalized Laguerre quadrature.
  359. Computes the sample points and weights for Gauss-generalized Laguerre
  360. quadrature. The sample points are the roots of the n-th degree generalized
  361. Laguerre polynomial, :math:`L^{\alpha}_n(x)`. These sample points and
  362. weights correctly integrate polynomials of degree :math:`2n - 1` or less
  363. over the interval :math:`[0, \infty]` with weight function
  364. :math:`f(x) = x^{\alpha} e^{-x}`.
  365. Parameters
  366. ----------
  367. n : int
  368. quadrature order
  369. alpha : float
  370. alpha must be > -1
  371. mu : bool, optional
  372. If True, return the sum of the weights, optional.
  373. Returns
  374. -------
  375. x : ndarray
  376. Sample points
  377. w : ndarray
  378. Weights
  379. mu : float
  380. Sum of the weights
  381. See Also
  382. --------
  383. scipy.integrate.quadrature
  384. scipy.integrate.fixed_quad
  385. """
  386. m = int(n)
  387. if n < 1 or n != m:
  388. raise ValueError("n must be a positive integer.")
  389. if alpha < -1:
  390. raise ValueError("alpha must be greater than -1.")
  391. mu0 = cephes.gamma(alpha + 1)
  392. if m == 1:
  393. x = np.array([alpha+1.0], 'd')
  394. w = np.array([mu0], 'd')
  395. if mu:
  396. return x, w, mu0
  397. else:
  398. return x, w
  399. an_func = lambda k: 2 * k + alpha + 1
  400. bn_func = lambda k: -np.sqrt(k * (k + alpha))
  401. f = lambda n, x: cephes.eval_genlaguerre(n, alpha, x)
  402. df = lambda n, x: (n*cephes.eval_genlaguerre(n, alpha, x)
  403. - (n + alpha)*cephes.eval_genlaguerre(n-1, alpha, x))/x
  404. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, False, mu)
  405. def genlaguerre(n, alpha, monic=False):
  406. r"""Generalized (associated) Laguerre polynomial.
  407. Defined to be the solution of
  408. .. math::
  409. x\frac{d^2}{dx^2}L_n^{(\alpha)}
  410. + (\alpha + 1 - x)\frac{d}{dx}L_n^{(\alpha)}
  411. + nL_n^{(\alpha)} = 0,
  412. where :math:`\alpha > -1`; :math:`L_n^{(\alpha)}` is a polynomial
  413. of degree :math:`n`.
  414. Parameters
  415. ----------
  416. n : int
  417. Degree of the polynomial.
  418. alpha : float
  419. Parameter, must be greater than -1.
  420. monic : bool, optional
  421. If `True`, scale the leading coefficient to be 1. Default is
  422. `False`.
  423. Returns
  424. -------
  425. L : orthopoly1d
  426. Generalized Laguerre polynomial.
  427. Notes
  428. -----
  429. For fixed :math:`\alpha`, the polynomials :math:`L_n^{(\alpha)}`
  430. are orthogonal over :math:`[0, \infty)` with weight function
  431. :math:`e^{-x}x^\alpha`.
  432. The Laguerre polynomials are the special case where :math:`\alpha
  433. = 0`.
  434. See Also
  435. --------
  436. laguerre : Laguerre polynomial.
  437. """
  438. if alpha <= -1:
  439. raise ValueError("alpha must be > -1")
  440. if n < 0:
  441. raise ValueError("n must be nonnegative.")
  442. if n == 0:
  443. n1 = n + 1
  444. else:
  445. n1 = n
  446. x, w, mu0 = roots_genlaguerre(n1, alpha, mu=True)
  447. wfunc = lambda x: exp(-x) * x**alpha
  448. if n == 0:
  449. x, w = [], []
  450. hn = _gam(n + alpha + 1) / _gam(n + 1)
  451. kn = (-1)**n / _gam(n + 1)
  452. p = orthopoly1d(x, w, hn, kn, wfunc, (0, inf), monic,
  453. lambda x: eval_genlaguerre(n, alpha, x))
  454. return p
  455. # Laguerre L_n(x)
  456. def roots_laguerre(n, mu=False):
  457. r"""Gauss-Laguerre quadrature.
  458. Computes the sample points and weights for Gauss-Laguerre quadrature.
  459. The sample points are the roots of the n-th degree Laguerre polynomial,
  460. :math:`L_n(x)`. These sample points and weights correctly integrate
  461. polynomials of degree :math:`2n - 1` or less over the interval
  462. :math:`[0, \infty]` with weight function :math:`f(x) = e^{-x}`.
  463. Parameters
  464. ----------
  465. n : int
  466. quadrature order
  467. mu : bool, optional
  468. If True, return the sum of the weights, optional.
  469. Returns
  470. -------
  471. x : ndarray
  472. Sample points
  473. w : ndarray
  474. Weights
  475. mu : float
  476. Sum of the weights
  477. See Also
  478. --------
  479. scipy.integrate.quadrature
  480. scipy.integrate.fixed_quad
  481. numpy.polynomial.laguerre.laggauss
  482. """
  483. return roots_genlaguerre(n, 0.0, mu=mu)
  484. def laguerre(n, monic=False):
  485. r"""Laguerre polynomial.
  486. Defined to be the solution of
  487. .. math::
  488. x\frac{d^2}{dx^2}L_n + (1 - x)\frac{d}{dx}L_n + nL_n = 0;
  489. :math:`L_n` is a polynomial of degree :math:`n`.
  490. Parameters
  491. ----------
  492. n : int
  493. Degree of the polynomial.
  494. monic : bool, optional
  495. If `True`, scale the leading coefficient to be 1. Default is
  496. `False`.
  497. Returns
  498. -------
  499. L : orthopoly1d
  500. Laguerre Polynomial.
  501. Notes
  502. -----
  503. The polynomials :math:`L_n` are orthogonal over :math:`[0,
  504. \infty)` with weight function :math:`e^{-x}`.
  505. """
  506. if n < 0:
  507. raise ValueError("n must be nonnegative.")
  508. if n == 0:
  509. n1 = n + 1
  510. else:
  511. n1 = n
  512. x, w, mu0 = roots_laguerre(n1, mu=True)
  513. if n == 0:
  514. x, w = [], []
  515. hn = 1.0
  516. kn = (-1)**n / _gam(n + 1)
  517. p = orthopoly1d(x, w, hn, kn, lambda x: exp(-x), (0, inf), monic,
  518. lambda x: eval_laguerre(n, x))
  519. return p
  520. # Hermite 1 H_n(x)
  521. def roots_hermite(n, mu=False):
  522. r"""Gauss-Hermite (physicst's) quadrature.
  523. Computes the sample points and weights for Gauss-Hermite quadrature.
  524. The sample points are the roots of the n-th degree Hermite polynomial,
  525. :math:`H_n(x)`. These sample points and weights correctly integrate
  526. polynomials of degree :math:`2n - 1` or less over the interval
  527. :math:`[-\infty, \infty]` with weight function :math:`f(x) = e^{-x^2}`.
  528. Parameters
  529. ----------
  530. n : int
  531. quadrature order
  532. mu : bool, optional
  533. If True, return the sum of the weights, optional.
  534. Returns
  535. -------
  536. x : ndarray
  537. Sample points
  538. w : ndarray
  539. Weights
  540. mu : float
  541. Sum of the weights
  542. Notes
  543. -----
  544. For small n up to 150 a modified version of the Golub-Welsch
  545. algorithm is used. Nodes are computed from the eigenvalue
  546. problem and improved by one step of a Newton iteration.
  547. The weights are computed from the well-known analytical formula.
  548. For n larger than 150 an optimal asymptotic algorithm is applied
  549. which computes nodes and weights in a numerically stable manner.
  550. The algorithm has linear runtime making computation for very
  551. large n (several thousand or more) feasible.
  552. See Also
  553. --------
  554. scipy.integrate.quadrature
  555. scipy.integrate.fixed_quad
  556. numpy.polynomial.hermite.hermgauss
  557. roots_hermitenorm
  558. References
  559. ----------
  560. .. [townsend.trogdon.olver-2014]
  561. Townsend, A. and Trogdon, T. and Olver, S. (2014)
  562. *Fast computation of Gauss quadrature nodes and
  563. weights on the whole real line*. :arXiv:`1410.5286`.
  564. .. [townsend.trogdon.olver-2015]
  565. Townsend, A. and Trogdon, T. and Olver, S. (2015)
  566. *Fast computation of Gauss quadrature nodes and
  567. weights on the whole real line*.
  568. IMA Journal of Numerical Analysis
  569. :doi:`10.1093/imanum/drv002`.
  570. """
  571. m = int(n)
  572. if n < 1 or n != m:
  573. raise ValueError("n must be a positive integer.")
  574. mu0 = np.sqrt(np.pi)
  575. if n <= 150:
  576. an_func = lambda k: 0.0*k
  577. bn_func = lambda k: np.sqrt(k/2.0)
  578. f = cephes.eval_hermite
  579. df = lambda n, x: 2.0 * n * cephes.eval_hermite(n-1, x)
  580. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  581. else:
  582. nodes, weights = _roots_hermite_asy(m)
  583. if mu:
  584. return nodes, weights, mu0
  585. else:
  586. return nodes, weights
  587. def _compute_tauk(n, k, maxit=5):
  588. """Helper function for Tricomi initial guesses
  589. For details, see formula 3.1 in lemma 3.1 in the
  590. original paper.
  591. Parameters
  592. ----------
  593. n : int
  594. Quadrature order
  595. k : ndarray of type int
  596. Index of roots :math:`\tau_k` to compute
  597. maxit : int
  598. Number of Newton maxit performed, the default
  599. value of 5 is sufficient.
  600. Returns
  601. -------
  602. tauk : ndarray
  603. Roots of equation 3.1
  604. See Also
  605. --------
  606. initial_nodes_a
  607. roots_hermite_asy
  608. """
  609. a = n % 2 - 0.5
  610. c = (4.0*floor(n/2.0) - 4.0*k + 3.0)*pi / (4.0*floor(n/2.0) + 2.0*a + 2.0)
  611. f = lambda x: x - sin(x) - c
  612. df = lambda x: 1.0 - cos(x)
  613. xi = 0.5*pi
  614. for i in range(maxit):
  615. xi = xi - f(xi)/df(xi)
  616. return xi
  617. def _initial_nodes_a(n, k):
  618. r"""Tricomi initial guesses
  619. Computes an initial approximation to the square of the `k`-th
  620. (positive) root :math:`x_k` of the Hermite polynomial :math:`H_n`
  621. of order :math:`n`. The formula is the one from lemma 3.1 in the
  622. original paper. The guesses are accurate except in the region
  623. near :math:`\sqrt{2n + 1}`.
  624. Parameters
  625. ----------
  626. n : int
  627. Quadrature order
  628. k : ndarray of type int
  629. Index of roots to compute
  630. Returns
  631. -------
  632. xksq : ndarray
  633. Square of the approximate roots
  634. See Also
  635. --------
  636. initial_nodes
  637. roots_hermite_asy
  638. """
  639. tauk = _compute_tauk(n, k)
  640. sigk = cos(0.5*tauk)**2
  641. a = n % 2 - 0.5
  642. nu = 4.0*floor(n/2.0) + 2.0*a + 2.0
  643. # Initial approximation of Hermite roots (square)
  644. xksq = nu*sigk - 1.0/(3.0*nu) * (5.0/(4.0*(1.0-sigk)**2) - 1.0/(1.0-sigk) - 0.25)
  645. return xksq
  646. def _initial_nodes_b(n, k):
  647. r"""Gatteschi initial guesses
  648. Computes an initial approximation to the square of the `k`-th
  649. (positive) root :math:`x_k` of the Hermite polynomial :math:`H_n`
  650. of order :math:`n`. The formula is the one from lemma 3.2 in the
  651. original paper. The guesses are accurate in the region just
  652. below :math:`\sqrt{2n + 1}`.
  653. Parameters
  654. ----------
  655. n : int
  656. Quadrature order
  657. k : ndarray of type int
  658. Index of roots to compute
  659. Returns
  660. -------
  661. xksq : ndarray
  662. Square of the approximate root
  663. See Also
  664. --------
  665. initial_nodes
  666. roots_hermite_asy
  667. """
  668. a = n % 2 - 0.5
  669. nu = 4.0*floor(n/2.0) + 2.0*a + 2.0
  670. # Airy roots by approximation
  671. ak = specfun.airyzo(k.max(), 1)[0][::-1]
  672. # Initial approximation of Hermite roots (square)
  673. xksq = (nu +
  674. 2.0**(2.0/3.0) * ak * nu**(1.0/3.0) +
  675. 1.0/5.0 * 2.0**(4.0/3.0) * ak**2 * nu**(-1.0/3.0) +
  676. (9.0/140.0 - 12.0/175.0 * ak**3) * nu**(-1.0) +
  677. (16.0/1575.0 * ak + 92.0/7875.0 * ak**4) * 2.0**(2.0/3.0) * nu**(-5.0/3.0) -
  678. (15152.0/3031875.0 * ak**5 + 1088.0/121275.0 * ak**2) * 2.0**(1.0/3.0) * nu**(-7.0/3.0))
  679. return xksq
  680. def _initial_nodes(n):
  681. """Initial guesses for the Hermite roots
  682. Computes an initial approximation to the non-negative
  683. roots :math:`x_k` of the Hermite polynomial :math:`H_n`
  684. of order :math:`n`. The Tricomi and Gatteschi initial
  685. guesses are used in the region where they are accurate.
  686. Parameters
  687. ----------
  688. n : int
  689. Quadrature order
  690. Returns
  691. -------
  692. xk : ndarray
  693. Approximate roots
  694. See Also
  695. --------
  696. roots_hermite_asy
  697. """
  698. # Turnover point
  699. # linear polynomial fit to error of 10, 25, 40, ..., 1000 point rules
  700. fit = 0.49082003*n - 4.37859653
  701. turnover = around(fit).astype(int)
  702. # Compute all approximations
  703. ia = arange(1, int(floor(n*0.5)+1))
  704. ib = ia[::-1]
  705. xasq = _initial_nodes_a(n, ia[:turnover+1])
  706. xbsq = _initial_nodes_b(n, ib[turnover+1:])
  707. # Combine
  708. iv = sqrt(hstack([xasq, xbsq]))
  709. # Central node is always zero
  710. if n % 2 == 1:
  711. iv = hstack([0.0, iv])
  712. return iv
  713. def _pbcf(n, theta):
  714. r"""Asymptotic series expansion of parabolic cylinder function
  715. The implementation is based on sections 3.2 and 3.3 from the
  716. original paper. Compared to the published version this code
  717. adds one more term to the asymptotic series. The detailed
  718. formulas can be found at [parabolic-asymptotics]_. The evaluation
  719. is done in a transformed variable :math:`\theta := \arccos(t)`
  720. where :math:`t := x / \mu` and :math:`\mu := \sqrt{2n + 1}`.
  721. Parameters
  722. ----------
  723. n : int
  724. Quadrature order
  725. theta : ndarray
  726. Transformed position variable
  727. Returns
  728. -------
  729. U : ndarray
  730. Value of the parabolic cylinder function :math:`U(a, \theta)`.
  731. Ud : ndarray
  732. Value of the derivative :math:`U^{\prime}(a, \theta)` of
  733. the parabolic cylinder function.
  734. See Also
  735. --------
  736. roots_hermite_asy
  737. References
  738. ----------
  739. .. [parabolic-asymptotics]
  740. https://dlmf.nist.gov/12.10#vii
  741. """
  742. st = sin(theta)
  743. ct = cos(theta)
  744. # https://dlmf.nist.gov/12.10#vii
  745. mu = 2.0*n + 1.0
  746. # https://dlmf.nist.gov/12.10#E23
  747. eta = 0.5*theta - 0.5*st*ct
  748. # https://dlmf.nist.gov/12.10#E39
  749. zeta = -(3.0*eta/2.0) ** (2.0/3.0)
  750. # https://dlmf.nist.gov/12.10#E40
  751. phi = (-zeta / st**2) ** (0.25)
  752. # Coefficients
  753. # https://dlmf.nist.gov/12.10#E43
  754. a0 = 1.0
  755. a1 = 0.10416666666666666667
  756. a2 = 0.08355034722222222222
  757. a3 = 0.12822657455632716049
  758. a4 = 0.29184902646414046425
  759. a5 = 0.88162726744375765242
  760. b0 = 1.0
  761. b1 = -0.14583333333333333333
  762. b2 = -0.09874131944444444444
  763. b3 = -0.14331205391589506173
  764. b4 = -0.31722720267841354810
  765. b5 = -0.94242914795712024914
  766. # Polynomials
  767. # https://dlmf.nist.gov/12.10#E9
  768. # https://dlmf.nist.gov/12.10#E10
  769. ctp = ct ** arange(16).reshape((-1,1))
  770. u0 = 1.0
  771. u1 = (1.0*ctp[3,:] - 6.0*ct) / 24.0
  772. u2 = (-9.0*ctp[4,:] + 249.0*ctp[2,:] + 145.0) / 1152.0
  773. u3 = (-4042.0*ctp[9,:] + 18189.0*ctp[7,:] - 28287.0*ctp[5,:] - 151995.0*ctp[3,:] - 259290.0*ct) / 414720.0
  774. u4 = (72756.0*ctp[10,:] - 321339.0*ctp[8,:] - 154982.0*ctp[6,:] + 50938215.0*ctp[4,:] + 122602962.0*ctp[2,:] + 12773113.0) / 39813120.0
  775. u5 = (82393456.0*ctp[15,:] - 617950920.0*ctp[13,:] + 1994971575.0*ctp[11,:] - 3630137104.0*ctp[9,:] + 4433574213.0*ctp[7,:]
  776. - 37370295816.0*ctp[5,:] - 119582875013.0*ctp[3,:] - 34009066266.0*ct) / 6688604160.0
  777. v0 = 1.0
  778. v1 = (1.0*ctp[3,:] + 6.0*ct) / 24.0
  779. v2 = (15.0*ctp[4,:] - 327.0*ctp[2,:] - 143.0) / 1152.0
  780. v3 = (-4042.0*ctp[9,:] + 18189.0*ctp[7,:] - 36387.0*ctp[5,:] + 238425.0*ctp[3,:] + 259290.0*ct) / 414720.0
  781. v4 = (-121260.0*ctp[10,:] + 551733.0*ctp[8,:] - 151958.0*ctp[6,:] - 57484425.0*ctp[4,:] - 132752238.0*ctp[2,:] - 12118727) / 39813120.0
  782. v5 = (82393456.0*ctp[15,:] - 617950920.0*ctp[13,:] + 2025529095.0*ctp[11,:] - 3750839308.0*ctp[9,:] + 3832454253.0*ctp[7,:]
  783. + 35213253348.0*ctp[5,:] + 130919230435.0*ctp[3,:] + 34009066266*ct) / 6688604160.0
  784. # Airy Evaluation (Bi and Bip unused)
  785. Ai, Aip, Bi, Bip = airy(mu**(4.0/6.0) * zeta)
  786. # Prefactor for U
  787. P = 2.0*sqrt(pi) * mu**(1.0/6.0) * phi
  788. # Terms for U
  789. # https://dlmf.nist.gov/12.10#E42
  790. phip = phi ** arange(6, 31, 6).reshape((-1,1))
  791. A0 = b0*u0
  792. A1 = (b2*u0 + phip[0,:]*b1*u1 + phip[1,:]*b0*u2) / zeta**3
  793. A2 = (b4*u0 + phip[0,:]*b3*u1 + phip[1,:]*b2*u2 + phip[2,:]*b1*u3 + phip[3,:]*b0*u4) / zeta**6
  794. B0 = -(a1*u0 + phip[0,:]*a0*u1) / zeta**2
  795. B1 = -(a3*u0 + phip[0,:]*a2*u1 + phip[1,:]*a1*u2 + phip[2,:]*a0*u3) / zeta**5
  796. B2 = -(a5*u0 + phip[0,:]*a4*u1 + phip[1,:]*a3*u2 + phip[2,:]*a2*u3 + phip[3,:]*a1*u4 + phip[4,:]*a0*u5) / zeta**8
  797. # U
  798. # https://dlmf.nist.gov/12.10#E35
  799. U = P * (Ai * (A0 + A1/mu**2.0 + A2/mu**4.0) +
  800. Aip * (B0 + B1/mu**2.0 + B2/mu**4.0) / mu**(8.0/6.0))
  801. # Prefactor for derivative of U
  802. Pd = sqrt(2.0*pi) * mu**(2.0/6.0) / phi
  803. # Terms for derivative of U
  804. # https://dlmf.nist.gov/12.10#E46
  805. C0 = -(b1*v0 + phip[0,:]*b0*v1) / zeta
  806. C1 = -(b3*v0 + phip[0,:]*b2*v1 + phip[1,:]*b1*v2 + phip[2,:]*b0*v3) / zeta**4
  807. C2 = -(b5*v0 + phip[0,:]*b4*v1 + phip[1,:]*b3*v2 + phip[2,:]*b2*v3 + phip[3,:]*b1*v4 + phip[4,:]*b0*v5) / zeta**7
  808. D0 = a0*v0
  809. D1 = (a2*v0 + phip[0,:]*a1*v1 + phip[1,:]*a0*v2) / zeta**3
  810. D2 = (a4*v0 + phip[0,:]*a3*v1 + phip[1,:]*a2*v2 + phip[2,:]*a1*v3 + phip[3,:]*a0*v4) / zeta**6
  811. # Derivative of U
  812. # https://dlmf.nist.gov/12.10#E36
  813. Ud = Pd * (Ai * (C0 + C1/mu**2.0 + C2/mu**4.0) / mu**(4.0/6.0) +
  814. Aip * (D0 + D1/mu**2.0 + D2/mu**4.0))
  815. return U, Ud
  816. def _newton(n, x_initial, maxit=5):
  817. """Newton iteration for polishing the asymptotic approximation
  818. to the zeros of the Hermite polynomials.
  819. Parameters
  820. ----------
  821. n : int
  822. Quadrature order
  823. x_initial : ndarray
  824. Initial guesses for the roots
  825. maxit : int
  826. Maximal number of Newton iterations.
  827. The default 5 is sufficient, usually
  828. only one or two steps are needed.
  829. Returns
  830. -------
  831. nodes : ndarray
  832. Quadrature nodes
  833. weights : ndarray
  834. Quadrature weights
  835. See Also
  836. --------
  837. roots_hermite_asy
  838. """
  839. # Variable transformation
  840. mu = sqrt(2.0*n + 1.0)
  841. t = x_initial / mu
  842. theta = arccos(t)
  843. # Newton iteration
  844. for i in range(maxit):
  845. u, ud = _pbcf(n, theta)
  846. dtheta = u / (sqrt(2.0) * mu * sin(theta) * ud)
  847. theta = theta + dtheta
  848. if max(abs(dtheta)) < 1e-14:
  849. break
  850. # Undo variable transformation
  851. x = mu * cos(theta)
  852. # Central node is always zero
  853. if n % 2 == 1:
  854. x[0] = 0.0
  855. # Compute weights
  856. w = exp(-x**2) / (2.0*ud**2)
  857. return x, w
  858. def _roots_hermite_asy(n):
  859. r"""Gauss-Hermite (physicst's) quadrature for large n.
  860. Computes the sample points and weights for Gauss-Hermite quadrature.
  861. The sample points are the roots of the n-th degree Hermite polynomial,
  862. :math:`H_n(x)`. These sample points and weights correctly integrate
  863. polynomials of degree :math:`2n - 1` or less over the interval
  864. :math:`[-\infty, \infty]` with weight function :math:`f(x) = e^{-x^2}`.
  865. This method relies on asymptotic expansions which work best for n > 150.
  866. The algorithm has linear runtime making computation for very large n
  867. feasible.
  868. Parameters
  869. ----------
  870. n : int
  871. quadrature order
  872. Returns
  873. -------
  874. nodes : ndarray
  875. Quadrature nodes
  876. weights : ndarray
  877. Quadrature weights
  878. See Also
  879. --------
  880. roots_hermite
  881. References
  882. ----------
  883. .. [townsend.trogdon.olver-2014]
  884. Townsend, A. and Trogdon, T. and Olver, S. (2014)
  885. *Fast computation of Gauss quadrature nodes and
  886. weights on the whole real line*. :arXiv:`1410.5286`.
  887. .. [townsend.trogdon.olver-2015]
  888. Townsend, A. and Trogdon, T. and Olver, S. (2015)
  889. *Fast computation of Gauss quadrature nodes and
  890. weights on the whole real line*.
  891. IMA Journal of Numerical Analysis
  892. :doi:`10.1093/imanum/drv002`.
  893. """
  894. iv = _initial_nodes(n)
  895. nodes, weights = _newton(n, iv)
  896. # Combine with negative parts
  897. if n % 2 == 0:
  898. nodes = hstack([-nodes[::-1], nodes])
  899. weights = hstack([weights[::-1], weights])
  900. else:
  901. nodes = hstack([-nodes[-1:0:-1], nodes])
  902. weights = hstack([weights[-1:0:-1], weights])
  903. # Scale weights
  904. weights *= sqrt(pi) / sum(weights)
  905. return nodes, weights
  906. def hermite(n, monic=False):
  907. r"""Physicist's Hermite polynomial.
  908. Defined by
  909. .. math::
  910. H_n(x) = (-1)^ne^{x^2}\frac{d^n}{dx^n}e^{-x^2};
  911. :math:`H_n` is a polynomial of degree :math:`n`.
  912. Parameters
  913. ----------
  914. n : int
  915. Degree of the polynomial.
  916. monic : bool, optional
  917. If `True`, scale the leading coefficient to be 1. Default is
  918. `False`.
  919. Returns
  920. -------
  921. H : orthopoly1d
  922. Hermite polynomial.
  923. Notes
  924. -----
  925. The polynomials :math:`H_n` are orthogonal over :math:`(-\infty,
  926. \infty)` with weight function :math:`e^{-x^2}`.
  927. """
  928. if n < 0:
  929. raise ValueError("n must be nonnegative.")
  930. if n == 0:
  931. n1 = n + 1
  932. else:
  933. n1 = n
  934. x, w, mu0 = roots_hermite(n1, mu=True)
  935. wfunc = lambda x: exp(-x * x)
  936. if n == 0:
  937. x, w = [], []
  938. hn = 2**n * _gam(n + 1) * sqrt(pi)
  939. kn = 2**n
  940. p = orthopoly1d(x, w, hn, kn, wfunc, (-inf, inf), monic,
  941. lambda x: eval_hermite(n, x))
  942. return p
  943. # Hermite 2 He_n(x)
  944. def roots_hermitenorm(n, mu=False):
  945. r"""Gauss-Hermite (statistician's) quadrature.
  946. Computes the sample points and weights for Gauss-Hermite quadrature.
  947. The sample points are the roots of the n-th degree Hermite polynomial,
  948. :math:`He_n(x)`. These sample points and weights correctly integrate
  949. polynomials of degree :math:`2n - 1` or less over the interval
  950. :math:`[-\infty, \infty]` with weight function :math:`f(x) = e^{-x^2/2}`.
  951. Parameters
  952. ----------
  953. n : int
  954. quadrature order
  955. mu : bool, optional
  956. If True, return the sum of the weights, optional.
  957. Returns
  958. -------
  959. x : ndarray
  960. Sample points
  961. w : ndarray
  962. Weights
  963. mu : float
  964. Sum of the weights
  965. Notes
  966. -----
  967. For small n up to 150 a modified version of the Golub-Welsch
  968. algorithm is used. Nodes are computed from the eigenvalue
  969. problem and improved by one step of a Newton iteration.
  970. The weights are computed from the well-known analytical formula.
  971. For n larger than 150 an optimal asymptotic algorithm is used
  972. which computes nodes and weights in a numerical stable manner.
  973. The algorithm has linear runtime making computation for very
  974. large n (several thousand or more) feasible.
  975. See Also
  976. --------
  977. scipy.integrate.quadrature
  978. scipy.integrate.fixed_quad
  979. numpy.polynomial.hermite_e.hermegauss
  980. """
  981. m = int(n)
  982. if n < 1 or n != m:
  983. raise ValueError("n must be a positive integer.")
  984. mu0 = np.sqrt(2.0*np.pi)
  985. if n <= 150:
  986. an_func = lambda k: 0.0*k
  987. bn_func = lambda k: np.sqrt(k)
  988. f = cephes.eval_hermitenorm
  989. df = lambda n, x: n * cephes.eval_hermitenorm(n-1, x)
  990. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  991. else:
  992. nodes, weights = _roots_hermite_asy(m)
  993. # Transform
  994. nodes *= sqrt(2)
  995. weights *= sqrt(2)
  996. if mu:
  997. return nodes, weights, mu0
  998. else:
  999. return nodes, weights
  1000. def hermitenorm(n, monic=False):
  1001. r"""Normalized (probabilist's) Hermite polynomial.
  1002. Defined by
  1003. .. math::
  1004. He_n(x) = (-1)^ne^{x^2/2}\frac{d^n}{dx^n}e^{-x^2/2};
  1005. :math:`He_n` is a polynomial of degree :math:`n`.
  1006. Parameters
  1007. ----------
  1008. n : int
  1009. Degree of the polynomial.
  1010. monic : bool, optional
  1011. If `True`, scale the leading coefficient to be 1. Default is
  1012. `False`.
  1013. Returns
  1014. -------
  1015. He : orthopoly1d
  1016. Hermite polynomial.
  1017. Notes
  1018. -----
  1019. The polynomials :math:`He_n` are orthogonal over :math:`(-\infty,
  1020. \infty)` with weight function :math:`e^{-x^2/2}`.
  1021. """
  1022. if n < 0:
  1023. raise ValueError("n must be nonnegative.")
  1024. if n == 0:
  1025. n1 = n + 1
  1026. else:
  1027. n1 = n
  1028. x, w, mu0 = roots_hermitenorm(n1, mu=True)
  1029. wfunc = lambda x: exp(-x * x / 2.0)
  1030. if n == 0:
  1031. x, w = [], []
  1032. hn = sqrt(2 * pi) * _gam(n + 1)
  1033. kn = 1.0
  1034. p = orthopoly1d(x, w, hn, kn, wfunc=wfunc, limits=(-inf, inf), monic=monic,
  1035. eval_func=lambda x: eval_hermitenorm(n, x))
  1036. return p
  1037. # The remainder of the polynomials can be derived from the ones above.
  1038. # Ultraspherical (Gegenbauer) C^(alpha)_n(x)
  1039. def roots_gegenbauer(n, alpha, mu=False):
  1040. r"""Gauss-Gegenbauer quadrature.
  1041. Computes the sample points and weights for Gauss-Gegenbauer quadrature.
  1042. The sample points are the roots of the n-th degree Gegenbauer polynomial,
  1043. :math:`C^{\alpha}_n(x)`. These sample points and weights correctly
  1044. integrate polynomials of degree :math:`2n - 1` or less over the interval
  1045. :math:`[-1, 1]` with weight function
  1046. :math:`f(x) = (1 - x^2)^{\alpha - 1/2}`.
  1047. Parameters
  1048. ----------
  1049. n : int
  1050. quadrature order
  1051. alpha : float
  1052. alpha must be > -0.5
  1053. mu : bool, optional
  1054. If True, return the sum of the weights, optional.
  1055. Returns
  1056. -------
  1057. x : ndarray
  1058. Sample points
  1059. w : ndarray
  1060. Weights
  1061. mu : float
  1062. Sum of the weights
  1063. See Also
  1064. --------
  1065. scipy.integrate.quadrature
  1066. scipy.integrate.fixed_quad
  1067. """
  1068. m = int(n)
  1069. if n < 1 or n != m:
  1070. raise ValueError("n must be a positive integer.")
  1071. if alpha < -0.5:
  1072. raise ValueError("alpha must be greater than -0.5.")
  1073. elif alpha == 0.0:
  1074. # C(n,0,x) == 0 uniformly, however, as alpha->0, C(n,alpha,x)->T(n,x)
  1075. # strictly, we should just error out here, since the roots are not
  1076. # really defined, but we used to return something useful, so let's
  1077. # keep doing so.
  1078. return roots_chebyt(n, mu)
  1079. mu0 = np.sqrt(np.pi) * cephes.gamma(alpha + 0.5) / cephes.gamma(alpha + 1)
  1080. an_func = lambda k: 0.0 * k
  1081. bn_func = lambda k: np.sqrt(k * (k + 2 * alpha - 1)
  1082. / (4 * (k + alpha) * (k + alpha - 1)))
  1083. f = lambda n, x: cephes.eval_gegenbauer(n, alpha, x)
  1084. df = lambda n, x: (-n*x*cephes.eval_gegenbauer(n, alpha, x)
  1085. + (n + 2*alpha - 1)*cephes.eval_gegenbauer(n-1, alpha, x))/(1-x**2)
  1086. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  1087. def gegenbauer(n, alpha, monic=False):
  1088. r"""Gegenbauer (ultraspherical) polynomial.
  1089. Defined to be the solution of
  1090. .. math::
  1091. (1 - x^2)\frac{d^2}{dx^2}C_n^{(\alpha)}
  1092. - (2\alpha + 1)x\frac{d}{dx}C_n^{(\alpha)}
  1093. + n(n + 2\alpha)C_n^{(\alpha)} = 0
  1094. for :math:`\alpha > -1/2`; :math:`C_n^{(\alpha)}` is a polynomial
  1095. of degree :math:`n`.
  1096. Parameters
  1097. ----------
  1098. n : int
  1099. Degree of the polynomial.
  1100. monic : bool, optional
  1101. If `True`, scale the leading coefficient to be 1. Default is
  1102. `False`.
  1103. Returns
  1104. -------
  1105. C : orthopoly1d
  1106. Gegenbauer polynomial.
  1107. Notes
  1108. -----
  1109. The polynomials :math:`C_n^{(\alpha)}` are orthogonal over
  1110. :math:`[-1,1]` with weight function :math:`(1 - x^2)^{(\alpha -
  1111. 1/2)}`.
  1112. """
  1113. base = jacobi(n, alpha - 0.5, alpha - 0.5, monic=monic)
  1114. if monic:
  1115. return base
  1116. # Abrahmowitz and Stegan 22.5.20
  1117. factor = (_gam(2*alpha + n) * _gam(alpha + 0.5) /
  1118. _gam(2*alpha) / _gam(alpha + 0.5 + n))
  1119. base._scale(factor)
  1120. base.__dict__['_eval_func'] = lambda x: eval_gegenbauer(float(n), alpha, x)
  1121. return base
  1122. # Chebyshev of the first kind: T_n(x) =
  1123. # n! sqrt(pi) / _gam(n+1./2)* P^(-1/2,-1/2)_n(x)
  1124. # Computed anew.
  1125. def roots_chebyt(n, mu=False):
  1126. r"""Gauss-Chebyshev (first kind) quadrature.
  1127. Computes the sample points and weights for Gauss-Chebyshev quadrature.
  1128. The sample points are the roots of the n-th degree Chebyshev polynomial of
  1129. the first kind, :math:`T_n(x)`. These sample points and weights correctly
  1130. integrate polynomials of degree :math:`2n - 1` or less over the interval
  1131. :math:`[-1, 1]` with weight function :math:`f(x) = 1/\sqrt{1 - x^2}`.
  1132. Parameters
  1133. ----------
  1134. n : int
  1135. quadrature order
  1136. mu : bool, optional
  1137. If True, return the sum of the weights, optional.
  1138. Returns
  1139. -------
  1140. x : ndarray
  1141. Sample points
  1142. w : ndarray
  1143. Weights
  1144. mu : float
  1145. Sum of the weights
  1146. See Also
  1147. --------
  1148. scipy.integrate.quadrature
  1149. scipy.integrate.fixed_quad
  1150. numpy.polynomial.chebyshev.chebgauss
  1151. """
  1152. m = int(n)
  1153. if n < 1 or n != m:
  1154. raise ValueError('n must be a positive integer.')
  1155. x = _ufuncs._sinpi(np.arange(-m + 1, m, 2) / (2*m))
  1156. w = np.full_like(x, pi/m)
  1157. if mu:
  1158. return x, w, pi
  1159. else:
  1160. return x, w
  1161. def chebyt(n, monic=False):
  1162. r"""Chebyshev polynomial of the first kind.
  1163. Defined to be the solution of
  1164. .. math::
  1165. (1 - x^2)\frac{d^2}{dx^2}T_n - x\frac{d}{dx}T_n + n^2T_n = 0;
  1166. :math:`T_n` is a polynomial of degree :math:`n`.
  1167. Parameters
  1168. ----------
  1169. n : int
  1170. Degree of the polynomial.
  1171. monic : bool, optional
  1172. If `True`, scale the leading coefficient to be 1. Default is
  1173. `False`.
  1174. Returns
  1175. -------
  1176. T : orthopoly1d
  1177. Chebyshev polynomial of the first kind.
  1178. Notes
  1179. -----
  1180. The polynomials :math:`T_n` are orthogonal over :math:`[-1, 1]`
  1181. with weight function :math:`(1 - x^2)^{-1/2}`.
  1182. See Also
  1183. --------
  1184. chebyu : Chebyshev polynomial of the second kind.
  1185. """
  1186. if n < 0:
  1187. raise ValueError("n must be nonnegative.")
  1188. wfunc = lambda x: 1.0 / sqrt(1 - x * x)
  1189. if n == 0:
  1190. return orthopoly1d([], [], pi, 1.0, wfunc, (-1, 1), monic,
  1191. lambda x: eval_chebyt(n, x))
  1192. n1 = n
  1193. x, w, mu = roots_chebyt(n1, mu=True)
  1194. hn = pi / 2
  1195. kn = 2**(n - 1)
  1196. p = orthopoly1d(x, w, hn, kn, wfunc, (-1, 1), monic,
  1197. lambda x: eval_chebyt(n, x))
  1198. return p
  1199. # Chebyshev of the second kind
  1200. # U_n(x) = (n+1)! sqrt(pi) / (2*_gam(n+3./2)) * P^(1/2,1/2)_n(x)
  1201. def roots_chebyu(n, mu=False):
  1202. r"""Gauss-Chebyshev (second kind) quadrature.
  1203. Computes the sample points and weights for Gauss-Chebyshev quadrature.
  1204. The sample points are the roots of the n-th degree Chebyshev polynomial of
  1205. the second kind, :math:`U_n(x)`. These sample points and weights correctly
  1206. integrate polynomials of degree :math:`2n - 1` or less over the interval
  1207. :math:`[-1, 1]` with weight function :math:`f(x) = \sqrt{1 - x^2}`.
  1208. Parameters
  1209. ----------
  1210. n : int
  1211. quadrature order
  1212. mu : bool, optional
  1213. If True, return the sum of the weights, optional.
  1214. Returns
  1215. -------
  1216. x : ndarray
  1217. Sample points
  1218. w : ndarray
  1219. Weights
  1220. mu : float
  1221. Sum of the weights
  1222. See Also
  1223. --------
  1224. scipy.integrate.quadrature
  1225. scipy.integrate.fixed_quad
  1226. """
  1227. m = int(n)
  1228. if n < 1 or n != m:
  1229. raise ValueError('n must be a positive integer.')
  1230. t = np.arange(m, 0, -1) * pi / (m + 1)
  1231. x = np.cos(t)
  1232. w = pi * np.sin(t)**2 / (m + 1)
  1233. if mu:
  1234. return x, w, pi / 2
  1235. else:
  1236. return x, w
  1237. def chebyu(n, monic=False):
  1238. r"""Chebyshev polynomial of the second kind.
  1239. Defined to be the solution of
  1240. .. math::
  1241. (1 - x^2)\frac{d^2}{dx^2}U_n - 3x\frac{d}{dx}U_n
  1242. + n(n + 2)U_n = 0;
  1243. :math:`U_n` is a polynomial of degree :math:`n`.
  1244. Parameters
  1245. ----------
  1246. n : int
  1247. Degree of the polynomial.
  1248. monic : bool, optional
  1249. If `True`, scale the leading coefficient to be 1. Default is
  1250. `False`.
  1251. Returns
  1252. -------
  1253. U : orthopoly1d
  1254. Chebyshev polynomial of the second kind.
  1255. Notes
  1256. -----
  1257. The polynomials :math:`U_n` are orthogonal over :math:`[-1, 1]`
  1258. with weight function :math:`(1 - x^2)^{1/2}`.
  1259. See Also
  1260. --------
  1261. chebyt : Chebyshev polynomial of the first kind.
  1262. """
  1263. base = jacobi(n, 0.5, 0.5, monic=monic)
  1264. if monic:
  1265. return base
  1266. factor = sqrt(pi) / 2.0 * _gam(n + 2) / _gam(n + 1.5)
  1267. base._scale(factor)
  1268. return base
  1269. # Chebyshev of the first kind C_n(x)
  1270. def roots_chebyc(n, mu=False):
  1271. r"""Gauss-Chebyshev (first kind) quadrature.
  1272. Computes the sample points and weights for Gauss-Chebyshev quadrature.
  1273. The sample points are the roots of the n-th degree Chebyshev polynomial of
  1274. the first kind, :math:`C_n(x)`. These sample points and weights correctly
  1275. integrate polynomials of degree :math:`2n - 1` or less over the interval
  1276. :math:`[-2, 2]` with weight function :math:`f(x) = 1/\sqrt{1 - (x/2)^2}`.
  1277. Parameters
  1278. ----------
  1279. n : int
  1280. quadrature order
  1281. mu : bool, optional
  1282. If True, return the sum of the weights, optional.
  1283. Returns
  1284. -------
  1285. x : ndarray
  1286. Sample points
  1287. w : ndarray
  1288. Weights
  1289. mu : float
  1290. Sum of the weights
  1291. See Also
  1292. --------
  1293. scipy.integrate.quadrature
  1294. scipy.integrate.fixed_quad
  1295. """
  1296. x, w, m = roots_chebyt(n, True)
  1297. x *= 2
  1298. w *= 2
  1299. m *= 2
  1300. if mu:
  1301. return x, w, m
  1302. else:
  1303. return x, w
  1304. def chebyc(n, monic=False):
  1305. r"""Chebyshev polynomial of the first kind on :math:`[-2, 2]`.
  1306. Defined as :math:`C_n(x) = 2T_n(x/2)`, where :math:`T_n` is the
  1307. nth Chebychev polynomial of the first kind.
  1308. Parameters
  1309. ----------
  1310. n : int
  1311. Degree of the polynomial.
  1312. monic : bool, optional
  1313. If `True`, scale the leading coefficient to be 1. Default is
  1314. `False`.
  1315. Returns
  1316. -------
  1317. C : orthopoly1d
  1318. Chebyshev polynomial of the first kind on :math:`[-2, 2]`.
  1319. Notes
  1320. -----
  1321. The polynomials :math:`C_n(x)` are orthogonal over :math:`[-2, 2]`
  1322. with weight function :math:`1/\sqrt{1 - (x/2)^2}`.
  1323. See Also
  1324. --------
  1325. chebyt : Chebyshev polynomial of the first kind.
  1326. References
  1327. ----------
  1328. .. [1] Abramowitz and Stegun, "Handbook of Mathematical Functions"
  1329. Section 22. National Bureau of Standards, 1972.
  1330. """
  1331. if n < 0:
  1332. raise ValueError("n must be nonnegative.")
  1333. if n == 0:
  1334. n1 = n + 1
  1335. else:
  1336. n1 = n
  1337. x, w, mu0 = roots_chebyc(n1, mu=True)
  1338. if n == 0:
  1339. x, w = [], []
  1340. hn = 4 * pi * ((n == 0) + 1)
  1341. kn = 1.0
  1342. p = orthopoly1d(x, w, hn, kn,
  1343. wfunc=lambda x: 1.0 / sqrt(1 - x * x / 4.0),
  1344. limits=(-2, 2), monic=monic)
  1345. if not monic:
  1346. p._scale(2.0 / p(2))
  1347. p.__dict__['_eval_func'] = lambda x: eval_chebyc(n, x)
  1348. return p
  1349. # Chebyshev of the second kind S_n(x)
  1350. def roots_chebys(n, mu=False):
  1351. r"""Gauss-Chebyshev (second kind) quadrature.
  1352. Computes the sample points and weights for Gauss-Chebyshev quadrature.
  1353. The sample points are the roots of the n-th degree Chebyshev polynomial of
  1354. the second kind, :math:`S_n(x)`. These sample points and weights correctly
  1355. integrate polynomials of degree :math:`2n - 1` or less over the interval
  1356. :math:`[-2, 2]` with weight function :math:`f(x) = \sqrt{1 - (x/2)^2}`.
  1357. Parameters
  1358. ----------
  1359. n : int
  1360. quadrature order
  1361. mu : bool, optional
  1362. If True, return the sum of the weights, optional.
  1363. Returns
  1364. -------
  1365. x : ndarray
  1366. Sample points
  1367. w : ndarray
  1368. Weights
  1369. mu : float
  1370. Sum of the weights
  1371. See Also
  1372. --------
  1373. scipy.integrate.quadrature
  1374. scipy.integrate.fixed_quad
  1375. """
  1376. x, w, m = roots_chebyu(n, True)
  1377. x *= 2
  1378. w *= 2
  1379. m *= 2
  1380. if mu:
  1381. return x, w, m
  1382. else:
  1383. return x, w
  1384. def chebys(n, monic=False):
  1385. r"""Chebyshev polynomial of the second kind on :math:`[-2, 2]`.
  1386. Defined as :math:`S_n(x) = U_n(x/2)` where :math:`U_n` is the
  1387. nth Chebychev polynomial of the second kind.
  1388. Parameters
  1389. ----------
  1390. n : int
  1391. Degree of the polynomial.
  1392. monic : bool, optional
  1393. If `True`, scale the leading coefficient to be 1. Default is
  1394. `False`.
  1395. Returns
  1396. -------
  1397. S : orthopoly1d
  1398. Chebyshev polynomial of the second kind on :math:`[-2, 2]`.
  1399. Notes
  1400. -----
  1401. The polynomials :math:`S_n(x)` are orthogonal over :math:`[-2, 2]`
  1402. with weight function :math:`\sqrt{1 - (x/2)}^2`.
  1403. See Also
  1404. --------
  1405. chebyu : Chebyshev polynomial of the second kind
  1406. References
  1407. ----------
  1408. .. [1] Abramowitz and Stegun, "Handbook of Mathematical Functions"
  1409. Section 22. National Bureau of Standards, 1972.
  1410. """
  1411. if n < 0:
  1412. raise ValueError("n must be nonnegative.")
  1413. if n == 0:
  1414. n1 = n + 1
  1415. else:
  1416. n1 = n
  1417. x, w, mu0 = roots_chebys(n1, mu=True)
  1418. if n == 0:
  1419. x, w = [], []
  1420. hn = pi
  1421. kn = 1.0
  1422. p = orthopoly1d(x, w, hn, kn,
  1423. wfunc=lambda x: sqrt(1 - x * x / 4.0),
  1424. limits=(-2, 2), monic=monic)
  1425. if not monic:
  1426. factor = (n + 1.0) / p(2)
  1427. p._scale(factor)
  1428. p.__dict__['_eval_func'] = lambda x: eval_chebys(n, x)
  1429. return p
  1430. # Shifted Chebyshev of the first kind T^*_n(x)
  1431. def roots_sh_chebyt(n, mu=False):
  1432. r"""Gauss-Chebyshev (first kind, shifted) quadrature.
  1433. Computes the sample points and weights for Gauss-Chebyshev quadrature.
  1434. The sample points are the roots of the n-th degree shifted Chebyshev
  1435. polynomial of the first kind, :math:`T_n(x)`. These sample points and
  1436. weights correctly integrate polynomials of degree :math:`2n - 1` or less
  1437. over the interval :math:`[0, 1]` with weight function
  1438. :math:`f(x) = 1/\sqrt{x - x^2}`.
  1439. Parameters
  1440. ----------
  1441. n : int
  1442. quadrature order
  1443. mu : bool, optional
  1444. If True, return the sum of the weights, optional.
  1445. Returns
  1446. -------
  1447. x : ndarray
  1448. Sample points
  1449. w : ndarray
  1450. Weights
  1451. mu : float
  1452. Sum of the weights
  1453. See Also
  1454. --------
  1455. scipy.integrate.quadrature
  1456. scipy.integrate.fixed_quad
  1457. """
  1458. xw = roots_chebyt(n, mu)
  1459. return ((xw[0] + 1) / 2,) + xw[1:]
  1460. def sh_chebyt(n, monic=False):
  1461. r"""Shifted Chebyshev polynomial of the first kind.
  1462. Defined as :math:`T^*_n(x) = T_n(2x - 1)` for :math:`T_n` the nth
  1463. Chebyshev polynomial of the first kind.
  1464. Parameters
  1465. ----------
  1466. n : int
  1467. Degree of the polynomial.
  1468. monic : bool, optional
  1469. If `True`, scale the leading coefficient to be 1. Default is
  1470. `False`.
  1471. Returns
  1472. -------
  1473. T : orthopoly1d
  1474. Shifted Chebyshev polynomial of the first kind.
  1475. Notes
  1476. -----
  1477. The polynomials :math:`T^*_n` are orthogonal over :math:`[0, 1]`
  1478. with weight function :math:`(x - x^2)^{-1/2}`.
  1479. """
  1480. base = sh_jacobi(n, 0.0, 0.5, monic=monic)
  1481. if monic:
  1482. return base
  1483. if n > 0:
  1484. factor = 4**n / 2.0
  1485. else:
  1486. factor = 1.0
  1487. base._scale(factor)
  1488. return base
  1489. # Shifted Chebyshev of the second kind U^*_n(x)
  1490. def roots_sh_chebyu(n, mu=False):
  1491. r"""Gauss-Chebyshev (second kind, shifted) quadrature.
  1492. Computes the sample points and weights for Gauss-Chebyshev quadrature.
  1493. The sample points are the roots of the n-th degree shifted Chebyshev
  1494. polynomial of the second kind, :math:`U_n(x)`. These sample points and
  1495. weights correctly integrate polynomials of degree :math:`2n - 1` or less
  1496. over the interval :math:`[0, 1]` with weight function
  1497. :math:`f(x) = \sqrt{x - x^2}`.
  1498. Parameters
  1499. ----------
  1500. n : int
  1501. quadrature order
  1502. mu : bool, optional
  1503. If True, return the sum of the weights, optional.
  1504. Returns
  1505. -------
  1506. x : ndarray
  1507. Sample points
  1508. w : ndarray
  1509. Weights
  1510. mu : float
  1511. Sum of the weights
  1512. See Also
  1513. --------
  1514. scipy.integrate.quadrature
  1515. scipy.integrate.fixed_quad
  1516. """
  1517. x, w, m = roots_chebyu(n, True)
  1518. x = (x + 1) / 2
  1519. m_us = cephes.beta(1.5, 1.5)
  1520. w *= m_us / m
  1521. if mu:
  1522. return x, w, m_us
  1523. else:
  1524. return x, w
  1525. def sh_chebyu(n, monic=False):
  1526. r"""Shifted Chebyshev polynomial of the second kind.
  1527. Defined as :math:`U^*_n(x) = U_n(2x - 1)` for :math:`U_n` the nth
  1528. Chebyshev polynomial of the second kind.
  1529. Parameters
  1530. ----------
  1531. n : int
  1532. Degree of the polynomial.
  1533. monic : bool, optional
  1534. If `True`, scale the leading coefficient to be 1. Default is
  1535. `False`.
  1536. Returns
  1537. -------
  1538. U : orthopoly1d
  1539. Shifted Chebyshev polynomial of the second kind.
  1540. Notes
  1541. -----
  1542. The polynomials :math:`U^*_n` are orthogonal over :math:`[0, 1]`
  1543. with weight function :math:`(x - x^2)^{1/2}`.
  1544. """
  1545. base = sh_jacobi(n, 2.0, 1.5, monic=monic)
  1546. if monic:
  1547. return base
  1548. factor = 4**n
  1549. base._scale(factor)
  1550. return base
  1551. # Legendre
  1552. def roots_legendre(n, mu=False):
  1553. r"""Gauss-Legendre quadrature.
  1554. Computes the sample points and weights for Gauss-Legendre quadrature.
  1555. The sample points are the roots of the n-th degree Legendre polynomial
  1556. :math:`P_n(x)`. These sample points and weights correctly integrate
  1557. polynomials of degree :math:`2n - 1` or less over the interval
  1558. :math:`[-1, 1]` with weight function :math:`f(x) = 1.0`.
  1559. Parameters
  1560. ----------
  1561. n : int
  1562. quadrature order
  1563. mu : bool, optional
  1564. If True, return the sum of the weights, optional.
  1565. Returns
  1566. -------
  1567. x : ndarray
  1568. Sample points
  1569. w : ndarray
  1570. Weights
  1571. mu : float
  1572. Sum of the weights
  1573. See Also
  1574. --------
  1575. scipy.integrate.quadrature
  1576. scipy.integrate.fixed_quad
  1577. numpy.polynomial.legendre.leggauss
  1578. """
  1579. m = int(n)
  1580. if n < 1 or n != m:
  1581. raise ValueError("n must be a positive integer.")
  1582. mu0 = 2.0
  1583. an_func = lambda k: 0.0 * k
  1584. bn_func = lambda k: k * np.sqrt(1.0 / (4 * k * k - 1))
  1585. f = cephes.eval_legendre
  1586. df = lambda n, x: (-n*x*cephes.eval_legendre(n, x)
  1587. + n*cephes.eval_legendre(n-1, x))/(1-x**2)
  1588. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  1589. def legendre(n, monic=False):
  1590. r"""Legendre polynomial.
  1591. Defined to be the solution of
  1592. .. math::
  1593. \frac{d}{dx}\left[(1 - x^2)\frac{d}{dx}P_n(x)\right]
  1594. + n(n + 1)P_n(x) = 0;
  1595. :math:`P_n(x)` is a polynomial of degree :math:`n`.
  1596. Parameters
  1597. ----------
  1598. n : int
  1599. Degree of the polynomial.
  1600. monic : bool, optional
  1601. If `True`, scale the leading coefficient to be 1. Default is
  1602. `False`.
  1603. Returns
  1604. -------
  1605. P : orthopoly1d
  1606. Legendre polynomial.
  1607. Notes
  1608. -----
  1609. The polynomials :math:`P_n` are orthogonal over :math:`[-1, 1]`
  1610. with weight function 1.
  1611. Examples
  1612. --------
  1613. Generate the 3rd-order Legendre polynomial 1/2*(5x^3 + 0x^2 - 3x + 0):
  1614. >>> from scipy.special import legendre
  1615. >>> legendre(3)
  1616. poly1d([ 2.5, 0. , -1.5, 0. ])
  1617. """
  1618. if n < 0:
  1619. raise ValueError("n must be nonnegative.")
  1620. if n == 0:
  1621. n1 = n + 1
  1622. else:
  1623. n1 = n
  1624. x, w, mu0 = roots_legendre(n1, mu=True)
  1625. if n == 0:
  1626. x, w = [], []
  1627. hn = 2.0 / (2 * n + 1)
  1628. kn = _gam(2 * n + 1) / _gam(n + 1)**2 / 2.0**n
  1629. p = orthopoly1d(x, w, hn, kn, wfunc=lambda x: 1.0, limits=(-1, 1),
  1630. monic=monic, eval_func=lambda x: eval_legendre(n, x))
  1631. return p
  1632. # Shifted Legendre P^*_n(x)
  1633. def roots_sh_legendre(n, mu=False):
  1634. r"""Gauss-Legendre (shifted) quadrature.
  1635. Computes the sample points and weights for Gauss-Legendre quadrature.
  1636. The sample points are the roots of the n-th degree shifted Legendre
  1637. polynomial :math:`P^*_n(x)`. These sample points and weights correctly
  1638. integrate polynomials of degree :math:`2n - 1` or less over the interval
  1639. :math:`[0, 1]` with weight function :math:`f(x) = 1.0`.
  1640. Parameters
  1641. ----------
  1642. n : int
  1643. quadrature order
  1644. mu : bool, optional
  1645. If True, return the sum of the weights, optional.
  1646. Returns
  1647. -------
  1648. x : ndarray
  1649. Sample points
  1650. w : ndarray
  1651. Weights
  1652. mu : float
  1653. Sum of the weights
  1654. See Also
  1655. --------
  1656. scipy.integrate.quadrature
  1657. scipy.integrate.fixed_quad
  1658. """
  1659. x, w = roots_legendre(n)
  1660. x = (x + 1) / 2
  1661. w /= 2
  1662. if mu:
  1663. return x, w, 1.0
  1664. else:
  1665. return x, w
  1666. def sh_legendre(n, monic=False):
  1667. r"""Shifted Legendre polynomial.
  1668. Defined as :math:`P^*_n(x) = P_n(2x - 1)` for :math:`P_n` the nth
  1669. Legendre polynomial.
  1670. Parameters
  1671. ----------
  1672. n : int
  1673. Degree of the polynomial.
  1674. monic : bool, optional
  1675. If `True`, scale the leading coefficient to be 1. Default is
  1676. `False`.
  1677. Returns
  1678. -------
  1679. P : orthopoly1d
  1680. Shifted Legendre polynomial.
  1681. Notes
  1682. -----
  1683. The polynomials :math:`P^*_n` are orthogonal over :math:`[0, 1]`
  1684. with weight function 1.
  1685. """
  1686. if n < 0:
  1687. raise ValueError("n must be nonnegative.")
  1688. wfunc = lambda x: 0.0 * x + 1.0
  1689. if n == 0:
  1690. return orthopoly1d([], [], 1.0, 1.0, wfunc, (0, 1), monic,
  1691. lambda x: eval_sh_legendre(n, x))
  1692. x, w, mu0 = roots_sh_legendre(n, mu=True)
  1693. hn = 1.0 / (2 * n + 1.0)
  1694. kn = _gam(2 * n + 1) / _gam(n + 1)**2
  1695. p = orthopoly1d(x, w, hn, kn, wfunc, limits=(0, 1), monic=monic,
  1696. eval_func=lambda x: eval_sh_legendre(n, x))
  1697. return p
  1698. # -----------------------------------------------------------------------------
  1699. # Code for backwards compatibility
  1700. # -----------------------------------------------------------------------------
  1701. # Import functions in case someone is still calling the orthogonal
  1702. # module directly. (They shouldn't be; it's not in the public API).
  1703. poch = cephes.poch
  1704. from ._ufuncs import (binom, eval_jacobi, eval_sh_jacobi, eval_gegenbauer,
  1705. eval_chebyt, eval_chebyu, eval_chebys, eval_chebyc,
  1706. eval_sh_chebyt, eval_sh_chebyu, eval_legendre,
  1707. eval_sh_legendre, eval_genlaguerre, eval_laguerre,
  1708. eval_hermite, eval_hermitenorm)
  1709. # Make the old root function names an alias for the new ones
  1710. _modattrs = globals()
  1711. for newfun, oldfun in _rootfuns_map.items():
  1712. _modattrs[oldfun] = _modattrs[newfun]
  1713. __all__.append(oldfun)