_orthogonal.py 72 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807180818091810181118121813181418151816181718181819182018211822182318241825182618271828182918301831183218331834183518361837183818391840184118421843184418451846184718481849185018511852185318541855185618571858185918601861186218631864186518661867186818691870187118721873187418751876187718781879188018811882188318841885188618871888188918901891189218931894189518961897189818991900190119021903190419051906190719081909191019111912191319141915191619171918191919201921192219231924192519261927192819291930193119321933193419351936193719381939194019411942194319441945194619471948194919501951195219531954195519561957195819591960196119621963196419651966196719681969197019711972197319741975197619771978197919801981198219831984198519861987198819891990199119921993199419951996199719981999200020012002200320042005200620072008200920102011201220132014201520162017201820192020202120222023202420252026202720282029203020312032203320342035203620372038203920402041204220432044204520462047204820492050205120522053205420552056205720582059206020612062206320642065206620672068206920702071207220732074207520762077207820792080208120822083208420852086208720882089209020912092209320942095209620972098209921002101210221032104210521062107210821092110211121122113211421152116211721182119212021212122212321242125212621272128212921302131213221332134213521362137213821392140214121422143214421452146214721482149215021512152215321542155215621572158215921602161216221632164216521662167216821692170217121722173217421752176217721782179218021812182218321842185218621872188218921902191219221932194219521962197219821992200220122022203220422052206220722082209221022112212221322142215221622172218221922202221222222232224222522262227222822292230223122322233223422352236223722382239224022412242224322442245224622472248224922502251225222532254225522562257225822592260226122622263226422652266226722682269227022712272227322742275227622772278227922802281228222832284228522862287228822892290229122922293229422952296229722982299230023012302230323042305230623072308230923102311231223132314231523162317231823192320232123222323232423252326232723282329233023312332233323342335233623372338233923402341234223432344234523462347234823492350235123522353235423552356235723582359236023612362236323642365236623672368236923702371237223732374237523762377237823792380238123822383238423852386238723882389239023912392239323942395239623972398239924002401240224032404240524062407240824092410241124122413241424152416241724182419242024212422242324242425242624272428242924302431243224332434243524362437243824392440244124422443244424452446244724482449245024512452245324542455245624572458245924602461246224632464246524662467246824692470247124722473247424752476247724782479248024812482248324842485248624872488248924902491249224932494249524962497249824992500250125022503250425052506250725082509251025112512251325142515251625172518251925202521252225232524252525262527252825292530253125322533253425352536253725382539254025412542254325442545254625472548254925502551255225532554255525562557
  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. # SciPy imports.
  55. import numpy as np
  56. from numpy import (exp, inf, pi, sqrt, floor, sin, cos, around,
  57. hstack, arccos, arange)
  58. from scipy import linalg
  59. from scipy.special import airy
  60. # Local imports.
  61. from . import _ufuncs
  62. _gam = _ufuncs.gamma
  63. # There is no .pyi file for _specfun
  64. from . import _specfun # type: ignore
  65. _polyfuns = ['legendre', 'chebyt', 'chebyu', 'chebyc', 'chebys',
  66. 'jacobi', 'laguerre', 'genlaguerre', 'hermite',
  67. 'hermitenorm', 'gegenbauer', 'sh_legendre', 'sh_chebyt',
  68. 'sh_chebyu', 'sh_jacobi']
  69. # Correspondence between new and old names of root functions
  70. _rootfuns_map = {'roots_legendre': 'p_roots',
  71. 'roots_chebyt': 't_roots',
  72. 'roots_chebyu': 'u_roots',
  73. 'roots_chebyc': 'c_roots',
  74. 'roots_chebys': 's_roots',
  75. 'roots_jacobi': 'j_roots',
  76. 'roots_laguerre': 'l_roots',
  77. 'roots_genlaguerre': 'la_roots',
  78. 'roots_hermite': 'h_roots',
  79. 'roots_hermitenorm': 'he_roots',
  80. 'roots_gegenbauer': 'cg_roots',
  81. 'roots_sh_legendre': 'ps_roots',
  82. 'roots_sh_chebyt': 'ts_roots',
  83. 'roots_sh_chebyu': 'us_roots',
  84. 'roots_sh_jacobi': 'js_roots'}
  85. __all__ = _polyfuns + list(_rootfuns_map.keys())
  86. class orthopoly1d(np.poly1d):
  87. def __init__(self, roots, weights=None, hn=1.0, kn=1.0, wfunc=None,
  88. limits=None, monic=False, eval_func=None):
  89. equiv_weights = [weights[k] / wfunc(roots[k]) for
  90. k in range(len(roots))]
  91. mu = sqrt(hn)
  92. if monic:
  93. evf = eval_func
  94. if evf:
  95. knn = kn
  96. eval_func = lambda x: evf(x) / knn
  97. mu = mu / abs(kn)
  98. kn = 1.0
  99. # compute coefficients from roots, then scale
  100. poly = np.poly1d(roots, r=True)
  101. np.poly1d.__init__(self, poly.coeffs * float(kn))
  102. self.weights = np.array(list(zip(roots, weights, equiv_weights)))
  103. self.weight_func = wfunc
  104. self.limits = limits
  105. self.normcoef = mu
  106. # Note: eval_func will be discarded on arithmetic
  107. self._eval_func = eval_func
  108. def __call__(self, v):
  109. if self._eval_func and not isinstance(v, np.poly1d):
  110. return self._eval_func(v)
  111. else:
  112. return np.poly1d.__call__(self, v)
  113. def _scale(self, p):
  114. if p == 1.0:
  115. return
  116. self._coeffs *= p
  117. evf = self._eval_func
  118. if evf:
  119. self._eval_func = lambda x: evf(x) * p
  120. self.normcoef *= p
  121. def _gen_roots_and_weights(n, mu0, an_func, bn_func, f, df, symmetrize, mu):
  122. """[x,w] = gen_roots_and_weights(n,an_func,sqrt_bn_func,mu)
  123. Returns the roots (x) of an nth order orthogonal polynomial,
  124. and weights (w) to use in appropriate Gaussian quadrature with that
  125. orthogonal polynomial.
  126. The polynomials have the recurrence relation
  127. P_n+1(x) = (x - A_n) P_n(x) - B_n P_n-1(x)
  128. an_func(n) should return A_n
  129. sqrt_bn_func(n) should return sqrt(B_n)
  130. mu ( = h_0 ) is the integral of the weight over the orthogonal
  131. interval
  132. """
  133. k = np.arange(n, dtype='d')
  134. c = np.zeros((2, n))
  135. c[0,1:] = bn_func(k[1:])
  136. c[1,:] = an_func(k)
  137. x = linalg.eigvals_banded(c, overwrite_a_band=True)
  138. # improve roots by one application of Newton's method
  139. y = f(n, x)
  140. dy = df(n, x)
  141. x -= y/dy
  142. # fm and dy may contain very large/small values, so we
  143. # log-normalize them to maintain precision in the product fm*dy
  144. fm = f(n-1, x)
  145. log_fm = np.log(np.abs(fm))
  146. log_dy = np.log(np.abs(dy))
  147. fm /= np.exp((log_fm.max() + log_fm.min()) / 2.)
  148. dy /= np.exp((log_dy.max() + log_dy.min()) / 2.)
  149. w = 1.0 / (fm * dy)
  150. if symmetrize:
  151. w = (w + w[::-1]) / 2
  152. x = (x - x[::-1]) / 2
  153. w *= mu0 / w.sum()
  154. if mu:
  155. return x, w, mu0
  156. else:
  157. return x, w
  158. # Jacobi Polynomials 1 P^(alpha,beta)_n(x)
  159. def roots_jacobi(n, alpha, beta, mu=False):
  160. r"""Gauss-Jacobi quadrature.
  161. Compute the sample points and weights for Gauss-Jacobi
  162. quadrature. The sample points are the roots of the nth degree
  163. Jacobi polynomial, :math:`P^{\alpha, \beta}_n(x)`. These sample
  164. points and weights correctly integrate polynomials of degree
  165. :math:`2n - 1` or less over the interval :math:`[-1, 1]` with
  166. weight function :math:`w(x) = (1 - x)^{\alpha} (1 +
  167. x)^{\beta}`. See 22.2.1 in [AS]_ for details.
  168. Parameters
  169. ----------
  170. n : int
  171. quadrature order
  172. alpha : float
  173. alpha must be > -1
  174. beta : float
  175. beta must be > -1
  176. mu : bool, optional
  177. If True, return the sum of the weights, optional.
  178. Returns
  179. -------
  180. x : ndarray
  181. Sample points
  182. w : ndarray
  183. Weights
  184. mu : float
  185. Sum of the weights
  186. See Also
  187. --------
  188. scipy.integrate.quadrature
  189. scipy.integrate.fixed_quad
  190. References
  191. ----------
  192. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  193. Handbook of Mathematical Functions with Formulas,
  194. Graphs, and Mathematical Tables. New York: Dover, 1972.
  195. """
  196. m = int(n)
  197. if n < 1 or n != m:
  198. raise ValueError("n must be a positive integer.")
  199. if alpha <= -1 or beta <= -1:
  200. raise ValueError("alpha and beta must be greater than -1.")
  201. if alpha == 0.0 and beta == 0.0:
  202. return roots_legendre(m, mu)
  203. if alpha == beta:
  204. return roots_gegenbauer(m, alpha+0.5, mu)
  205. if (alpha + beta) <= 1000:
  206. mu0 = 2.0**(alpha+beta+1) * _ufuncs.beta(alpha+1, beta+1)
  207. else:
  208. # Avoid overflows in pow and beta for very large parameters
  209. mu0 = np.exp((alpha + beta + 1) * np.log(2.0)
  210. + _ufuncs.betaln(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: _ufuncs.eval_jacobi(n, a, b, x)
  221. df = lambda n, x: (0.5 * (n + a + b + 1)
  222. * _ufuncs.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. References
  255. ----------
  256. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  257. Handbook of Mathematical Functions with Formulas,
  258. Graphs, and Mathematical Tables. New York: Dover, 1972.
  259. Examples
  260. --------
  261. The Jacobi polynomials satisfy the recurrence relation:
  262. .. math::
  263. P_n^{(\alpha, \beta-1)}(x) - P_n^{(\alpha-1, \beta)}(x)
  264. = P_{n-1}^{(\alpha, \beta)}(x)
  265. This can be verified, for example, for :math:`\alpha = \beta = 2`
  266. and :math:`n = 1` over the interval :math:`[-1, 1]`:
  267. >>> import numpy as np
  268. >>> from scipy.special import jacobi
  269. >>> x = np.arange(-1.0, 1.0, 0.01)
  270. >>> np.allclose(jacobi(0, 2, 2)(x),
  271. ... jacobi(1, 2, 1)(x) - jacobi(1, 1, 2)(x))
  272. True
  273. Plot of the Jacobi polynomial :math:`P_5^{(\alpha, -0.5)}` for
  274. different values of :math:`\alpha`:
  275. >>> import matplotlib.pyplot as plt
  276. >>> x = np.arange(-1.0, 1.0, 0.01)
  277. >>> fig, ax = plt.subplots()
  278. >>> ax.set_ylim(-2.0, 2.0)
  279. >>> ax.set_title(r'Jacobi polynomials $P_5^{(\alpha, -0.5)}$')
  280. >>> for alpha in np.arange(0, 4, 1):
  281. ... ax.plot(x, jacobi(5, alpha, -0.5)(x), label=rf'$\alpha={alpha}$')
  282. >>> plt.legend(loc='best')
  283. >>> plt.show()
  284. """
  285. if n < 0:
  286. raise ValueError("n must be nonnegative.")
  287. wfunc = lambda x: (1 - x)**alpha * (1 + x)**beta
  288. if n == 0:
  289. return orthopoly1d([], [], 1.0, 1.0, wfunc, (-1, 1), monic,
  290. eval_func=np.ones_like)
  291. x, w, mu = roots_jacobi(n, alpha, beta, mu=True)
  292. ab1 = alpha + beta + 1.0
  293. hn = 2**ab1 / (2 * n + ab1) * _gam(n + alpha + 1)
  294. hn *= _gam(n + beta + 1.0) / _gam(n + 1) / _gam(n + ab1)
  295. kn = _gam(2 * n + ab1) / 2.0**n / _gam(n + 1) / _gam(n + ab1)
  296. # here kn = coefficient on x^n term
  297. p = orthopoly1d(x, w, hn, kn, wfunc, (-1, 1), monic,
  298. lambda x: _ufuncs.eval_jacobi(n, alpha, beta, x))
  299. return p
  300. # Jacobi Polynomials shifted G_n(p,q,x)
  301. def roots_sh_jacobi(n, p1, q1, mu=False):
  302. """Gauss-Jacobi (shifted) quadrature.
  303. Compute the sample points and weights for Gauss-Jacobi (shifted)
  304. quadrature. The sample points are the roots of the nth degree
  305. shifted Jacobi polynomial, :math:`G^{p,q}_n(x)`. These sample
  306. points and weights correctly integrate polynomials of degree
  307. :math:`2n - 1` or less over the interval :math:`[0, 1]` with
  308. weight function :math:`w(x) = (1 - x)^{p-q} x^{q-1}`. See 22.2.2
  309. in [AS]_ for details.
  310. Parameters
  311. ----------
  312. n : int
  313. quadrature order
  314. p1 : float
  315. (p1 - q1) must be > -1
  316. q1 : float
  317. q1 must be > 0
  318. mu : bool, optional
  319. If True, return the sum of the weights, optional.
  320. Returns
  321. -------
  322. x : ndarray
  323. Sample points
  324. w : ndarray
  325. Weights
  326. mu : float
  327. Sum of the weights
  328. See Also
  329. --------
  330. scipy.integrate.quadrature
  331. scipy.integrate.fixed_quad
  332. References
  333. ----------
  334. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  335. Handbook of Mathematical Functions with Formulas,
  336. Graphs, and Mathematical Tables. New York: Dover, 1972.
  337. """
  338. if (p1-q1) <= -1 or q1 <= 0:
  339. raise ValueError("(p - q) must be greater than -1, and q must be greater than 0.")
  340. x, w, m = roots_jacobi(n, p1-q1, q1-1, True)
  341. x = (x + 1) / 2
  342. scale = 2.0**p1
  343. w /= scale
  344. m /= scale
  345. if mu:
  346. return x, w, m
  347. else:
  348. return x, w
  349. def sh_jacobi(n, p, q, monic=False):
  350. r"""Shifted Jacobi polynomial.
  351. Defined by
  352. .. math::
  353. G_n^{(p, q)}(x)
  354. = \binom{2n + p - 1}{n}^{-1}P_n^{(p - q, q - 1)}(2x - 1),
  355. where :math:`P_n^{(\cdot, \cdot)}` is the nth Jacobi polynomial.
  356. Parameters
  357. ----------
  358. n : int
  359. Degree of the polynomial.
  360. p : float
  361. Parameter, must have :math:`p > q - 1`.
  362. q : float
  363. Parameter, must be greater than 0.
  364. monic : bool, optional
  365. If `True`, scale the leading coefficient to be 1. Default is
  366. `False`.
  367. Returns
  368. -------
  369. G : orthopoly1d
  370. Shifted Jacobi polynomial.
  371. Notes
  372. -----
  373. For fixed :math:`p, q`, the polynomials :math:`G_n^{(p, q)}` are
  374. orthogonal over :math:`[0, 1]` with weight function :math:`(1 -
  375. x)^{p - q}x^{q - 1}`.
  376. """
  377. if n < 0:
  378. raise ValueError("n must be nonnegative.")
  379. wfunc = lambda x: (1.0 - x)**(p - q) * (x)**(q - 1.)
  380. if n == 0:
  381. return orthopoly1d([], [], 1.0, 1.0, wfunc, (-1, 1), monic,
  382. eval_func=np.ones_like)
  383. n1 = n
  384. x, w = roots_sh_jacobi(n1, p, q)
  385. hn = _gam(n + 1) * _gam(n + q) * _gam(n + p) * _gam(n + p - q + 1)
  386. hn /= (2 * n + p) * (_gam(2 * n + p)**2)
  387. # kn = 1.0 in standard form so monic is redundant. Kept for compatibility.
  388. kn = 1.0
  389. pp = orthopoly1d(x, w, hn, kn, wfunc=wfunc, limits=(0, 1), monic=monic,
  390. eval_func=lambda x: _ufuncs.eval_sh_jacobi(n, p, q, x))
  391. return pp
  392. # Generalized Laguerre L^(alpha)_n(x)
  393. def roots_genlaguerre(n, alpha, mu=False):
  394. r"""Gauss-generalized Laguerre quadrature.
  395. Compute the sample points and weights for Gauss-generalized
  396. Laguerre quadrature. The sample points are the roots of the nth
  397. degree generalized Laguerre polynomial, :math:`L^{\alpha}_n(x)`.
  398. These sample points and weights correctly integrate polynomials of
  399. degree :math:`2n - 1` or less over the interval :math:`[0,
  400. \infty]` with weight function :math:`w(x) = x^{\alpha}
  401. e^{-x}`. See 22.3.9 in [AS]_ for details.
  402. Parameters
  403. ----------
  404. n : int
  405. quadrature order
  406. alpha : float
  407. alpha must be > -1
  408. mu : bool, optional
  409. If True, return the sum of the weights, optional.
  410. Returns
  411. -------
  412. x : ndarray
  413. Sample points
  414. w : ndarray
  415. Weights
  416. mu : float
  417. Sum of the weights
  418. See Also
  419. --------
  420. scipy.integrate.quadrature
  421. scipy.integrate.fixed_quad
  422. References
  423. ----------
  424. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  425. Handbook of Mathematical Functions with Formulas,
  426. Graphs, and Mathematical Tables. New York: Dover, 1972.
  427. """
  428. m = int(n)
  429. if n < 1 or n != m:
  430. raise ValueError("n must be a positive integer.")
  431. if alpha < -1:
  432. raise ValueError("alpha must be greater than -1.")
  433. mu0 = _ufuncs.gamma(alpha + 1)
  434. if m == 1:
  435. x = np.array([alpha+1.0], 'd')
  436. w = np.array([mu0], 'd')
  437. if mu:
  438. return x, w, mu0
  439. else:
  440. return x, w
  441. an_func = lambda k: 2 * k + alpha + 1
  442. bn_func = lambda k: -np.sqrt(k * (k + alpha))
  443. f = lambda n, x: _ufuncs.eval_genlaguerre(n, alpha, x)
  444. df = lambda n, x: (n*_ufuncs.eval_genlaguerre(n, alpha, x)
  445. - (n + alpha)*_ufuncs.eval_genlaguerre(n-1, alpha, x))/x
  446. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, False, mu)
  447. def genlaguerre(n, alpha, monic=False):
  448. r"""Generalized (associated) Laguerre polynomial.
  449. Defined to be the solution of
  450. .. math::
  451. x\frac{d^2}{dx^2}L_n^{(\alpha)}
  452. + (\alpha + 1 - x)\frac{d}{dx}L_n^{(\alpha)}
  453. + nL_n^{(\alpha)} = 0,
  454. where :math:`\alpha > -1`; :math:`L_n^{(\alpha)}` is a polynomial
  455. of degree :math:`n`.
  456. Parameters
  457. ----------
  458. n : int
  459. Degree of the polynomial.
  460. alpha : float
  461. Parameter, must be greater than -1.
  462. monic : bool, optional
  463. If `True`, scale the leading coefficient to be 1. Default is
  464. `False`.
  465. Returns
  466. -------
  467. L : orthopoly1d
  468. Generalized Laguerre polynomial.
  469. Notes
  470. -----
  471. For fixed :math:`\alpha`, the polynomials :math:`L_n^{(\alpha)}`
  472. are orthogonal over :math:`[0, \infty)` with weight function
  473. :math:`e^{-x}x^\alpha`.
  474. The Laguerre polynomials are the special case where :math:`\alpha
  475. = 0`.
  476. See Also
  477. --------
  478. laguerre : Laguerre polynomial.
  479. hyp1f1 : confluent hypergeometric function
  480. References
  481. ----------
  482. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  483. Handbook of Mathematical Functions with Formulas,
  484. Graphs, and Mathematical Tables. New York: Dover, 1972.
  485. Examples
  486. --------
  487. The generalized Laguerre polynomials are closely related to the confluent
  488. hypergeometric function :math:`{}_1F_1`:
  489. .. math::
  490. L_n^{(\alpha)} = \binom{n + \alpha}{n} {}_1F_1(-n, \alpha +1, x)
  491. This can be verified, for example, for :math:`n = \alpha = 3` over the
  492. interval :math:`[-1, 1]`:
  493. >>> import numpy as np
  494. >>> from scipy.special import binom
  495. >>> from scipy.special import genlaguerre
  496. >>> from scipy.special import hyp1f1
  497. >>> x = np.arange(-1.0, 1.0, 0.01)
  498. >>> np.allclose(genlaguerre(3, 3)(x), binom(6, 3) * hyp1f1(-3, 4, x))
  499. True
  500. This is the plot of the generalized Laguerre polynomials
  501. :math:`L_3^{(\alpha)}` for some values of :math:`\alpha`:
  502. >>> import matplotlib.pyplot as plt
  503. >>> x = np.arange(-4.0, 12.0, 0.01)
  504. >>> fig, ax = plt.subplots()
  505. >>> ax.set_ylim(-5.0, 10.0)
  506. >>> ax.set_title(r'Generalized Laguerre polynomials $L_3^{\alpha}$')
  507. >>> for alpha in np.arange(0, 5):
  508. ... ax.plot(x, genlaguerre(3, alpha)(x), label=rf'$L_3^{(alpha)}$')
  509. >>> plt.legend(loc='best')
  510. >>> plt.show()
  511. """
  512. if alpha <= -1:
  513. raise ValueError("alpha must be > -1")
  514. if n < 0:
  515. raise ValueError("n must be nonnegative.")
  516. if n == 0:
  517. n1 = n + 1
  518. else:
  519. n1 = n
  520. x, w = roots_genlaguerre(n1, alpha)
  521. wfunc = lambda x: exp(-x) * x**alpha
  522. if n == 0:
  523. x, w = [], []
  524. hn = _gam(n + alpha + 1) / _gam(n + 1)
  525. kn = (-1)**n / _gam(n + 1)
  526. p = orthopoly1d(x, w, hn, kn, wfunc, (0, inf), monic,
  527. lambda x: _ufuncs.eval_genlaguerre(n, alpha, x))
  528. return p
  529. # Laguerre L_n(x)
  530. def roots_laguerre(n, mu=False):
  531. r"""Gauss-Laguerre quadrature.
  532. Compute the sample points and weights for Gauss-Laguerre
  533. quadrature. The sample points are the roots of the nth degree
  534. Laguerre polynomial, :math:`L_n(x)`. These sample points and
  535. weights correctly integrate polynomials of degree :math:`2n - 1`
  536. or less over the interval :math:`[0, \infty]` with weight function
  537. :math:`w(x) = e^{-x}`. See 22.2.13 in [AS]_ for details.
  538. Parameters
  539. ----------
  540. n : int
  541. quadrature order
  542. mu : bool, optional
  543. If True, return the sum of the weights, optional.
  544. Returns
  545. -------
  546. x : ndarray
  547. Sample points
  548. w : ndarray
  549. Weights
  550. mu : float
  551. Sum of the weights
  552. See Also
  553. --------
  554. scipy.integrate.quadrature
  555. scipy.integrate.fixed_quad
  556. numpy.polynomial.laguerre.laggauss
  557. References
  558. ----------
  559. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  560. Handbook of Mathematical Functions with Formulas,
  561. Graphs, and Mathematical Tables. New York: Dover, 1972.
  562. """
  563. return roots_genlaguerre(n, 0.0, mu=mu)
  564. def laguerre(n, monic=False):
  565. r"""Laguerre polynomial.
  566. Defined to be the solution of
  567. .. math::
  568. x\frac{d^2}{dx^2}L_n + (1 - x)\frac{d}{dx}L_n + nL_n = 0;
  569. :math:`L_n` is a polynomial of degree :math:`n`.
  570. Parameters
  571. ----------
  572. n : int
  573. Degree of the polynomial.
  574. monic : bool, optional
  575. If `True`, scale the leading coefficient to be 1. Default is
  576. `False`.
  577. Returns
  578. -------
  579. L : orthopoly1d
  580. Laguerre Polynomial.
  581. Notes
  582. -----
  583. The polynomials :math:`L_n` are orthogonal over :math:`[0,
  584. \infty)` with weight function :math:`e^{-x}`.
  585. See Also
  586. --------
  587. genlaguerre : Generalized (associated) Laguerre polynomial.
  588. References
  589. ----------
  590. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  591. Handbook of Mathematical Functions with Formulas,
  592. Graphs, and Mathematical Tables. New York: Dover, 1972.
  593. Examples
  594. --------
  595. The Laguerre polynomials :math:`L_n` are the special case
  596. :math:`\alpha = 0` of the generalized Laguerre polynomials
  597. :math:`L_n^{(\alpha)}`.
  598. Let's verify it on the interval :math:`[-1, 1]`:
  599. >>> import numpy as np
  600. >>> from scipy.special import genlaguerre
  601. >>> from scipy.special import laguerre
  602. >>> x = np.arange(-1.0, 1.0, 0.01)
  603. >>> np.allclose(genlaguerre(3, 0)(x), laguerre(3)(x))
  604. True
  605. The polynomials :math:`L_n` also satisfy the recurrence relation:
  606. .. math::
  607. (n + 1)L_{n+1}(x) = (2n +1 -x)L_n(x) - nL_{n-1}(x)
  608. This can be easily checked on :math:`[0, 1]` for :math:`n = 3`:
  609. >>> x = np.arange(0.0, 1.0, 0.01)
  610. >>> np.allclose(4 * laguerre(4)(x),
  611. ... (7 - x) * laguerre(3)(x) - 3 * laguerre(2)(x))
  612. True
  613. This is the plot of the first few Laguerre polynomials :math:`L_n`:
  614. >>> import matplotlib.pyplot as plt
  615. >>> x = np.arange(-1.0, 5.0, 0.01)
  616. >>> fig, ax = plt.subplots()
  617. >>> ax.set_ylim(-5.0, 5.0)
  618. >>> ax.set_title(r'Laguerre polynomials $L_n$')
  619. >>> for n in np.arange(0, 5):
  620. ... ax.plot(x, laguerre(n)(x), label=rf'$L_{n}$')
  621. >>> plt.legend(loc='best')
  622. >>> plt.show()
  623. """
  624. if n < 0:
  625. raise ValueError("n must be nonnegative.")
  626. if n == 0:
  627. n1 = n + 1
  628. else:
  629. n1 = n
  630. x, w = roots_laguerre(n1)
  631. if n == 0:
  632. x, w = [], []
  633. hn = 1.0
  634. kn = (-1)**n / _gam(n + 1)
  635. p = orthopoly1d(x, w, hn, kn, lambda x: exp(-x), (0, inf), monic,
  636. lambda x: _ufuncs.eval_laguerre(n, x))
  637. return p
  638. # Hermite 1 H_n(x)
  639. def roots_hermite(n, mu=False):
  640. r"""Gauss-Hermite (physicist's) quadrature.
  641. Compute the sample points and weights for Gauss-Hermite
  642. quadrature. The sample points are the roots of the nth degree
  643. Hermite polynomial, :math:`H_n(x)`. These sample points and
  644. weights correctly integrate polynomials of degree :math:`2n - 1`
  645. or less over the interval :math:`[-\infty, \infty]` with weight
  646. function :math:`w(x) = e^{-x^2}`. See 22.2.14 in [AS]_ for
  647. details.
  648. Parameters
  649. ----------
  650. n : int
  651. quadrature order
  652. mu : bool, optional
  653. If True, return the sum of the weights, optional.
  654. Returns
  655. -------
  656. x : ndarray
  657. Sample points
  658. w : ndarray
  659. Weights
  660. mu : float
  661. Sum of the weights
  662. Notes
  663. -----
  664. For small n up to 150 a modified version of the Golub-Welsch
  665. algorithm is used. Nodes are computed from the eigenvalue
  666. problem and improved by one step of a Newton iteration.
  667. The weights are computed from the well-known analytical formula.
  668. For n larger than 150 an optimal asymptotic algorithm is applied
  669. which computes nodes and weights in a numerically stable manner.
  670. The algorithm has linear runtime making computation for very
  671. large n (several thousand or more) feasible.
  672. See Also
  673. --------
  674. scipy.integrate.quadrature
  675. scipy.integrate.fixed_quad
  676. numpy.polynomial.hermite.hermgauss
  677. roots_hermitenorm
  678. References
  679. ----------
  680. .. [townsend.trogdon.olver-2014]
  681. Townsend, A. and Trogdon, T. and Olver, S. (2014)
  682. *Fast computation of Gauss quadrature nodes and
  683. weights on the whole real line*. :arXiv:`1410.5286`.
  684. .. [townsend.trogdon.olver-2015]
  685. Townsend, A. and Trogdon, T. and Olver, S. (2015)
  686. *Fast computation of Gauss quadrature nodes and
  687. weights on the whole real line*.
  688. IMA Journal of Numerical Analysis
  689. :doi:`10.1093/imanum/drv002`.
  690. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  691. Handbook of Mathematical Functions with Formulas,
  692. Graphs, and Mathematical Tables. New York: Dover, 1972.
  693. """
  694. m = int(n)
  695. if n < 1 or n != m:
  696. raise ValueError("n must be a positive integer.")
  697. mu0 = np.sqrt(np.pi)
  698. if n <= 150:
  699. an_func = lambda k: 0.0*k
  700. bn_func = lambda k: np.sqrt(k/2.0)
  701. f = _ufuncs.eval_hermite
  702. df = lambda n, x: 2.0 * n * _ufuncs.eval_hermite(n-1, x)
  703. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  704. else:
  705. nodes, weights = _roots_hermite_asy(m)
  706. if mu:
  707. return nodes, weights, mu0
  708. else:
  709. return nodes, weights
  710. def _compute_tauk(n, k, maxit=5):
  711. """Helper function for Tricomi initial guesses
  712. For details, see formula 3.1 in lemma 3.1 in the
  713. original paper.
  714. Parameters
  715. ----------
  716. n : int
  717. Quadrature order
  718. k : ndarray of type int
  719. Index of roots :math:`\tau_k` to compute
  720. maxit : int
  721. Number of Newton maxit performed, the default
  722. value of 5 is sufficient.
  723. Returns
  724. -------
  725. tauk : ndarray
  726. Roots of equation 3.1
  727. See Also
  728. --------
  729. initial_nodes_a
  730. roots_hermite_asy
  731. """
  732. a = n % 2 - 0.5
  733. c = (4.0*floor(n/2.0) - 4.0*k + 3.0)*pi / (4.0*floor(n/2.0) + 2.0*a + 2.0)
  734. f = lambda x: x - sin(x) - c
  735. df = lambda x: 1.0 - cos(x)
  736. xi = 0.5*pi
  737. for i in range(maxit):
  738. xi = xi - f(xi)/df(xi)
  739. return xi
  740. def _initial_nodes_a(n, k):
  741. r"""Tricomi initial guesses
  742. Computes an initial approximation to the square of the `k`-th
  743. (positive) root :math:`x_k` of the Hermite polynomial :math:`H_n`
  744. of order :math:`n`. The formula is the one from lemma 3.1 in the
  745. original paper. The guesses are accurate except in the region
  746. near :math:`\sqrt{2n + 1}`.
  747. Parameters
  748. ----------
  749. n : int
  750. Quadrature order
  751. k : ndarray of type int
  752. Index of roots to compute
  753. Returns
  754. -------
  755. xksq : ndarray
  756. Square of the approximate roots
  757. See Also
  758. --------
  759. initial_nodes
  760. roots_hermite_asy
  761. """
  762. tauk = _compute_tauk(n, k)
  763. sigk = cos(0.5*tauk)**2
  764. a = n % 2 - 0.5
  765. nu = 4.0*floor(n/2.0) + 2.0*a + 2.0
  766. # Initial approximation of Hermite roots (square)
  767. xksq = nu*sigk - 1.0/(3.0*nu) * (5.0/(4.0*(1.0-sigk)**2) - 1.0/(1.0-sigk) - 0.25)
  768. return xksq
  769. def _initial_nodes_b(n, k):
  770. r"""Gatteschi initial guesses
  771. Computes an initial approximation to the square of the kth
  772. (positive) root :math:`x_k` of the Hermite polynomial :math:`H_n`
  773. of order :math:`n`. The formula is the one from lemma 3.2 in the
  774. original paper. The guesses are accurate in the region just
  775. below :math:`\sqrt{2n + 1}`.
  776. Parameters
  777. ----------
  778. n : int
  779. Quadrature order
  780. k : ndarray of type int
  781. Index of roots to compute
  782. Returns
  783. -------
  784. xksq : ndarray
  785. Square of the approximate root
  786. See Also
  787. --------
  788. initial_nodes
  789. roots_hermite_asy
  790. """
  791. a = n % 2 - 0.5
  792. nu = 4.0*floor(n/2.0) + 2.0*a + 2.0
  793. # Airy roots by approximation
  794. ak = _specfun.airyzo(k.max(), 1)[0][::-1]
  795. # Initial approximation of Hermite roots (square)
  796. xksq = (nu +
  797. 2.0**(2.0/3.0) * ak * nu**(1.0/3.0) +
  798. 1.0/5.0 * 2.0**(4.0/3.0) * ak**2 * nu**(-1.0/3.0) +
  799. (9.0/140.0 - 12.0/175.0 * ak**3) * nu**(-1.0) +
  800. (16.0/1575.0 * ak + 92.0/7875.0 * ak**4) * 2.0**(2.0/3.0) * nu**(-5.0/3.0) -
  801. (15152.0/3031875.0 * ak**5 + 1088.0/121275.0 * ak**2) * 2.0**(1.0/3.0) * nu**(-7.0/3.0))
  802. return xksq
  803. def _initial_nodes(n):
  804. """Initial guesses for the Hermite roots
  805. Computes an initial approximation to the non-negative
  806. roots :math:`x_k` of the Hermite polynomial :math:`H_n`
  807. of order :math:`n`. The Tricomi and Gatteschi initial
  808. guesses are used in the region where they are accurate.
  809. Parameters
  810. ----------
  811. n : int
  812. Quadrature order
  813. Returns
  814. -------
  815. xk : ndarray
  816. Approximate roots
  817. See Also
  818. --------
  819. roots_hermite_asy
  820. """
  821. # Turnover point
  822. # linear polynomial fit to error of 10, 25, 40, ..., 1000 point rules
  823. fit = 0.49082003*n - 4.37859653
  824. turnover = around(fit).astype(int)
  825. # Compute all approximations
  826. ia = arange(1, int(floor(n*0.5)+1))
  827. ib = ia[::-1]
  828. xasq = _initial_nodes_a(n, ia[:turnover+1])
  829. xbsq = _initial_nodes_b(n, ib[turnover+1:])
  830. # Combine
  831. iv = sqrt(hstack([xasq, xbsq]))
  832. # Central node is always zero
  833. if n % 2 == 1:
  834. iv = hstack([0.0, iv])
  835. return iv
  836. def _pbcf(n, theta):
  837. r"""Asymptotic series expansion of parabolic cylinder function
  838. The implementation is based on sections 3.2 and 3.3 from the
  839. original paper. Compared to the published version this code
  840. adds one more term to the asymptotic series. The detailed
  841. formulas can be found at [parabolic-asymptotics]_. The evaluation
  842. is done in a transformed variable :math:`\theta := \arccos(t)`
  843. where :math:`t := x / \mu` and :math:`\mu := \sqrt{2n + 1}`.
  844. Parameters
  845. ----------
  846. n : int
  847. Quadrature order
  848. theta : ndarray
  849. Transformed position variable
  850. Returns
  851. -------
  852. U : ndarray
  853. Value of the parabolic cylinder function :math:`U(a, \theta)`.
  854. Ud : ndarray
  855. Value of the derivative :math:`U^{\prime}(a, \theta)` of
  856. the parabolic cylinder function.
  857. See Also
  858. --------
  859. roots_hermite_asy
  860. References
  861. ----------
  862. .. [parabolic-asymptotics]
  863. https://dlmf.nist.gov/12.10#vii
  864. """
  865. st = sin(theta)
  866. ct = cos(theta)
  867. # https://dlmf.nist.gov/12.10#vii
  868. mu = 2.0*n + 1.0
  869. # https://dlmf.nist.gov/12.10#E23
  870. eta = 0.5*theta - 0.5*st*ct
  871. # https://dlmf.nist.gov/12.10#E39
  872. zeta = -(3.0*eta/2.0) ** (2.0/3.0)
  873. # https://dlmf.nist.gov/12.10#E40
  874. phi = (-zeta / st**2) ** (0.25)
  875. # Coefficients
  876. # https://dlmf.nist.gov/12.10#E43
  877. a0 = 1.0
  878. a1 = 0.10416666666666666667
  879. a2 = 0.08355034722222222222
  880. a3 = 0.12822657455632716049
  881. a4 = 0.29184902646414046425
  882. a5 = 0.88162726744375765242
  883. b0 = 1.0
  884. b1 = -0.14583333333333333333
  885. b2 = -0.09874131944444444444
  886. b3 = -0.14331205391589506173
  887. b4 = -0.31722720267841354810
  888. b5 = -0.94242914795712024914
  889. # Polynomials
  890. # https://dlmf.nist.gov/12.10#E9
  891. # https://dlmf.nist.gov/12.10#E10
  892. ctp = ct ** arange(16).reshape((-1,1))
  893. u0 = 1.0
  894. u1 = (1.0*ctp[3,:] - 6.0*ct) / 24.0
  895. u2 = (-9.0*ctp[4,:] + 249.0*ctp[2,:] + 145.0) / 1152.0
  896. u3 = (-4042.0*ctp[9,:] + 18189.0*ctp[7,:] - 28287.0*ctp[5,:] - 151995.0*ctp[3,:] - 259290.0*ct) / 414720.0
  897. 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
  898. u5 = (82393456.0*ctp[15,:] - 617950920.0*ctp[13,:] + 1994971575.0*ctp[11,:] - 3630137104.0*ctp[9,:] + 4433574213.0*ctp[7,:]
  899. - 37370295816.0*ctp[5,:] - 119582875013.0*ctp[3,:] - 34009066266.0*ct) / 6688604160.0
  900. v0 = 1.0
  901. v1 = (1.0*ctp[3,:] + 6.0*ct) / 24.0
  902. v2 = (15.0*ctp[4,:] - 327.0*ctp[2,:] - 143.0) / 1152.0
  903. v3 = (-4042.0*ctp[9,:] + 18189.0*ctp[7,:] - 36387.0*ctp[5,:] + 238425.0*ctp[3,:] + 259290.0*ct) / 414720.0
  904. v4 = (-121260.0*ctp[10,:] + 551733.0*ctp[8,:] - 151958.0*ctp[6,:] - 57484425.0*ctp[4,:] - 132752238.0*ctp[2,:] - 12118727) / 39813120.0
  905. v5 = (82393456.0*ctp[15,:] - 617950920.0*ctp[13,:] + 2025529095.0*ctp[11,:] - 3750839308.0*ctp[9,:] + 3832454253.0*ctp[7,:]
  906. + 35213253348.0*ctp[5,:] + 130919230435.0*ctp[3,:] + 34009066266*ct) / 6688604160.0
  907. # Airy Evaluation (Bi and Bip unused)
  908. Ai, Aip, Bi, Bip = airy(mu**(4.0/6.0) * zeta)
  909. # Prefactor for U
  910. P = 2.0*sqrt(pi) * mu**(1.0/6.0) * phi
  911. # Terms for U
  912. # https://dlmf.nist.gov/12.10#E42
  913. phip = phi ** arange(6, 31, 6).reshape((-1,1))
  914. A0 = b0*u0
  915. A1 = (b2*u0 + phip[0,:]*b1*u1 + phip[1,:]*b0*u2) / zeta**3
  916. A2 = (b4*u0 + phip[0,:]*b3*u1 + phip[1,:]*b2*u2 + phip[2,:]*b1*u3 + phip[3,:]*b0*u4) / zeta**6
  917. B0 = -(a1*u0 + phip[0,:]*a0*u1) / zeta**2
  918. B1 = -(a3*u0 + phip[0,:]*a2*u1 + phip[1,:]*a1*u2 + phip[2,:]*a0*u3) / zeta**5
  919. 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
  920. # U
  921. # https://dlmf.nist.gov/12.10#E35
  922. U = P * (Ai * (A0 + A1/mu**2.0 + A2/mu**4.0) +
  923. Aip * (B0 + B1/mu**2.0 + B2/mu**4.0) / mu**(8.0/6.0))
  924. # Prefactor for derivative of U
  925. Pd = sqrt(2.0*pi) * mu**(2.0/6.0) / phi
  926. # Terms for derivative of U
  927. # https://dlmf.nist.gov/12.10#E46
  928. C0 = -(b1*v0 + phip[0,:]*b0*v1) / zeta
  929. C1 = -(b3*v0 + phip[0,:]*b2*v1 + phip[1,:]*b1*v2 + phip[2,:]*b0*v3) / zeta**4
  930. 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
  931. D0 = a0*v0
  932. D1 = (a2*v0 + phip[0,:]*a1*v1 + phip[1,:]*a0*v2) / zeta**3
  933. D2 = (a4*v0 + phip[0,:]*a3*v1 + phip[1,:]*a2*v2 + phip[2,:]*a1*v3 + phip[3,:]*a0*v4) / zeta**6
  934. # Derivative of U
  935. # https://dlmf.nist.gov/12.10#E36
  936. Ud = Pd * (Ai * (C0 + C1/mu**2.0 + C2/mu**4.0) / mu**(4.0/6.0) +
  937. Aip * (D0 + D1/mu**2.0 + D2/mu**4.0))
  938. return U, Ud
  939. def _newton(n, x_initial, maxit=5):
  940. """Newton iteration for polishing the asymptotic approximation
  941. to the zeros of the Hermite polynomials.
  942. Parameters
  943. ----------
  944. n : int
  945. Quadrature order
  946. x_initial : ndarray
  947. Initial guesses for the roots
  948. maxit : int
  949. Maximal number of Newton iterations.
  950. The default 5 is sufficient, usually
  951. only one or two steps are needed.
  952. Returns
  953. -------
  954. nodes : ndarray
  955. Quadrature nodes
  956. weights : ndarray
  957. Quadrature weights
  958. See Also
  959. --------
  960. roots_hermite_asy
  961. """
  962. # Variable transformation
  963. mu = sqrt(2.0*n + 1.0)
  964. t = x_initial / mu
  965. theta = arccos(t)
  966. # Newton iteration
  967. for i in range(maxit):
  968. u, ud = _pbcf(n, theta)
  969. dtheta = u / (sqrt(2.0) * mu * sin(theta) * ud)
  970. theta = theta + dtheta
  971. if max(abs(dtheta)) < 1e-14:
  972. break
  973. # Undo variable transformation
  974. x = mu * cos(theta)
  975. # Central node is always zero
  976. if n % 2 == 1:
  977. x[0] = 0.0
  978. # Compute weights
  979. w = exp(-x**2) / (2.0*ud**2)
  980. return x, w
  981. def _roots_hermite_asy(n):
  982. r"""Gauss-Hermite (physicist's) quadrature for large n.
  983. Computes the sample points and weights for Gauss-Hermite quadrature.
  984. The sample points are the roots of the nth degree Hermite polynomial,
  985. :math:`H_n(x)`. These sample points and weights correctly integrate
  986. polynomials of degree :math:`2n - 1` or less over the interval
  987. :math:`[-\infty, \infty]` with weight function :math:`f(x) = e^{-x^2}`.
  988. This method relies on asymptotic expansions which work best for n > 150.
  989. The algorithm has linear runtime making computation for very large n
  990. feasible.
  991. Parameters
  992. ----------
  993. n : int
  994. quadrature order
  995. Returns
  996. -------
  997. nodes : ndarray
  998. Quadrature nodes
  999. weights : ndarray
  1000. Quadrature weights
  1001. See Also
  1002. --------
  1003. roots_hermite
  1004. References
  1005. ----------
  1006. .. [townsend.trogdon.olver-2014]
  1007. Townsend, A. and Trogdon, T. and Olver, S. (2014)
  1008. *Fast computation of Gauss quadrature nodes and
  1009. weights on the whole real line*. :arXiv:`1410.5286`.
  1010. .. [townsend.trogdon.olver-2015]
  1011. Townsend, A. and Trogdon, T. and Olver, S. (2015)
  1012. *Fast computation of Gauss quadrature nodes and
  1013. weights on the whole real line*.
  1014. IMA Journal of Numerical Analysis
  1015. :doi:`10.1093/imanum/drv002`.
  1016. """
  1017. iv = _initial_nodes(n)
  1018. nodes, weights = _newton(n, iv)
  1019. # Combine with negative parts
  1020. if n % 2 == 0:
  1021. nodes = hstack([-nodes[::-1], nodes])
  1022. weights = hstack([weights[::-1], weights])
  1023. else:
  1024. nodes = hstack([-nodes[-1:0:-1], nodes])
  1025. weights = hstack([weights[-1:0:-1], weights])
  1026. # Scale weights
  1027. weights *= sqrt(pi) / sum(weights)
  1028. return nodes, weights
  1029. def hermite(n, monic=False):
  1030. r"""Physicist's Hermite polynomial.
  1031. Defined by
  1032. .. math::
  1033. H_n(x) = (-1)^ne^{x^2}\frac{d^n}{dx^n}e^{-x^2};
  1034. :math:`H_n` is a polynomial of degree :math:`n`.
  1035. Parameters
  1036. ----------
  1037. n : int
  1038. Degree of the polynomial.
  1039. monic : bool, optional
  1040. If `True`, scale the leading coefficient to be 1. Default is
  1041. `False`.
  1042. Returns
  1043. -------
  1044. H : orthopoly1d
  1045. Hermite polynomial.
  1046. Notes
  1047. -----
  1048. The polynomials :math:`H_n` are orthogonal over :math:`(-\infty,
  1049. \infty)` with weight function :math:`e^{-x^2}`.
  1050. Examples
  1051. --------
  1052. >>> from scipy import special
  1053. >>> import matplotlib.pyplot as plt
  1054. >>> import numpy as np
  1055. >>> p_monic = special.hermite(3, monic=True)
  1056. >>> p_monic
  1057. poly1d([ 1. , 0. , -1.5, 0. ])
  1058. >>> p_monic(1)
  1059. -0.49999999999999983
  1060. >>> x = np.linspace(-3, 3, 400)
  1061. >>> y = p_monic(x)
  1062. >>> plt.plot(x, y)
  1063. >>> plt.title("Monic Hermite polynomial of degree 3")
  1064. >>> plt.xlabel("x")
  1065. >>> plt.ylabel("H_3(x)")
  1066. >>> plt.show()
  1067. """
  1068. if n < 0:
  1069. raise ValueError("n must be nonnegative.")
  1070. if n == 0:
  1071. n1 = n + 1
  1072. else:
  1073. n1 = n
  1074. x, w = roots_hermite(n1)
  1075. wfunc = lambda x: exp(-x * x)
  1076. if n == 0:
  1077. x, w = [], []
  1078. hn = 2**n * _gam(n + 1) * sqrt(pi)
  1079. kn = 2**n
  1080. p = orthopoly1d(x, w, hn, kn, wfunc, (-inf, inf), monic,
  1081. lambda x: _ufuncs.eval_hermite(n, x))
  1082. return p
  1083. # Hermite 2 He_n(x)
  1084. def roots_hermitenorm(n, mu=False):
  1085. r"""Gauss-Hermite (statistician's) quadrature.
  1086. Compute the sample points and weights for Gauss-Hermite
  1087. quadrature. The sample points are the roots of the nth degree
  1088. Hermite polynomial, :math:`He_n(x)`. These sample points and
  1089. weights correctly integrate polynomials of degree :math:`2n - 1`
  1090. or less over the interval :math:`[-\infty, \infty]` with weight
  1091. function :math:`w(x) = e^{-x^2/2}`. See 22.2.15 in [AS]_ for more
  1092. details.
  1093. Parameters
  1094. ----------
  1095. n : int
  1096. quadrature order
  1097. mu : bool, optional
  1098. If True, return the sum of the weights, optional.
  1099. Returns
  1100. -------
  1101. x : ndarray
  1102. Sample points
  1103. w : ndarray
  1104. Weights
  1105. mu : float
  1106. Sum of the weights
  1107. Notes
  1108. -----
  1109. For small n up to 150 a modified version of the Golub-Welsch
  1110. algorithm is used. Nodes are computed from the eigenvalue
  1111. problem and improved by one step of a Newton iteration.
  1112. The weights are computed from the well-known analytical formula.
  1113. For n larger than 150 an optimal asymptotic algorithm is used
  1114. which computes nodes and weights in a numerical stable manner.
  1115. The algorithm has linear runtime making computation for very
  1116. large n (several thousand or more) feasible.
  1117. See Also
  1118. --------
  1119. scipy.integrate.quadrature
  1120. scipy.integrate.fixed_quad
  1121. numpy.polynomial.hermite_e.hermegauss
  1122. References
  1123. ----------
  1124. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1125. Handbook of Mathematical Functions with Formulas,
  1126. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1127. """
  1128. m = int(n)
  1129. if n < 1 or n != m:
  1130. raise ValueError("n must be a positive integer.")
  1131. mu0 = np.sqrt(2.0*np.pi)
  1132. if n <= 150:
  1133. an_func = lambda k: 0.0*k
  1134. bn_func = lambda k: np.sqrt(k)
  1135. f = _ufuncs.eval_hermitenorm
  1136. df = lambda n, x: n * _ufuncs.eval_hermitenorm(n-1, x)
  1137. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  1138. else:
  1139. nodes, weights = _roots_hermite_asy(m)
  1140. # Transform
  1141. nodes *= sqrt(2)
  1142. weights *= sqrt(2)
  1143. if mu:
  1144. return nodes, weights, mu0
  1145. else:
  1146. return nodes, weights
  1147. def hermitenorm(n, monic=False):
  1148. r"""Normalized (probabilist's) Hermite polynomial.
  1149. Defined by
  1150. .. math::
  1151. He_n(x) = (-1)^ne^{x^2/2}\frac{d^n}{dx^n}e^{-x^2/2};
  1152. :math:`He_n` is a polynomial of degree :math:`n`.
  1153. Parameters
  1154. ----------
  1155. n : int
  1156. Degree of the polynomial.
  1157. monic : bool, optional
  1158. If `True`, scale the leading coefficient to be 1. Default is
  1159. `False`.
  1160. Returns
  1161. -------
  1162. He : orthopoly1d
  1163. Hermite polynomial.
  1164. Notes
  1165. -----
  1166. The polynomials :math:`He_n` are orthogonal over :math:`(-\infty,
  1167. \infty)` with weight function :math:`e^{-x^2/2}`.
  1168. """
  1169. if n < 0:
  1170. raise ValueError("n must be nonnegative.")
  1171. if n == 0:
  1172. n1 = n + 1
  1173. else:
  1174. n1 = n
  1175. x, w = roots_hermitenorm(n1)
  1176. wfunc = lambda x: exp(-x * x / 2.0)
  1177. if n == 0:
  1178. x, w = [], []
  1179. hn = sqrt(2 * pi) * _gam(n + 1)
  1180. kn = 1.0
  1181. p = orthopoly1d(x, w, hn, kn, wfunc=wfunc, limits=(-inf, inf), monic=monic,
  1182. eval_func=lambda x: _ufuncs.eval_hermitenorm(n, x))
  1183. return p
  1184. # The remainder of the polynomials can be derived from the ones above.
  1185. # Ultraspherical (Gegenbauer) C^(alpha)_n(x)
  1186. def roots_gegenbauer(n, alpha, mu=False):
  1187. r"""Gauss-Gegenbauer quadrature.
  1188. Compute the sample points and weights for Gauss-Gegenbauer
  1189. quadrature. The sample points are the roots of the nth degree
  1190. Gegenbauer polynomial, :math:`C^{\alpha}_n(x)`. These sample
  1191. points and weights correctly integrate polynomials of degree
  1192. :math:`2n - 1` or less over the interval :math:`[-1, 1]` with
  1193. weight function :math:`w(x) = (1 - x^2)^{\alpha - 1/2}`. See
  1194. 22.2.3 in [AS]_ for more details.
  1195. Parameters
  1196. ----------
  1197. n : int
  1198. quadrature order
  1199. alpha : float
  1200. alpha must be > -0.5
  1201. mu : bool, optional
  1202. If True, return the sum of the weights, optional.
  1203. Returns
  1204. -------
  1205. x : ndarray
  1206. Sample points
  1207. w : ndarray
  1208. Weights
  1209. mu : float
  1210. Sum of the weights
  1211. See Also
  1212. --------
  1213. scipy.integrate.quadrature
  1214. scipy.integrate.fixed_quad
  1215. References
  1216. ----------
  1217. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1218. Handbook of Mathematical Functions with Formulas,
  1219. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1220. """
  1221. m = int(n)
  1222. if n < 1 or n != m:
  1223. raise ValueError("n must be a positive integer.")
  1224. if alpha < -0.5:
  1225. raise ValueError("alpha must be greater than -0.5.")
  1226. elif alpha == 0.0:
  1227. # C(n,0,x) == 0 uniformly, however, as alpha->0, C(n,alpha,x)->T(n,x)
  1228. # strictly, we should just error out here, since the roots are not
  1229. # really defined, but we used to return something useful, so let's
  1230. # keep doing so.
  1231. return roots_chebyt(n, mu)
  1232. if alpha <= 170:
  1233. mu0 = (np.sqrt(np.pi) * _ufuncs.gamma(alpha + 0.5)) \
  1234. / _ufuncs.gamma(alpha + 1)
  1235. else:
  1236. # For large alpha we use a Taylor series expansion around inf,
  1237. # expressed as a 6th order polynomial of a^-1 and using Horner's
  1238. # method to minimize computation and maximize precision
  1239. inv_alpha = 1. / alpha
  1240. coeffs = np.array([0.000207186, -0.00152206, -0.000640869,
  1241. 0.00488281, 0.0078125, -0.125, 1.])
  1242. mu0 = coeffs[0]
  1243. for term in range(1, len(coeffs)):
  1244. mu0 = mu0 * inv_alpha + coeffs[term]
  1245. mu0 = mu0 * np.sqrt(np.pi / alpha)
  1246. an_func = lambda k: 0.0 * k
  1247. bn_func = lambda k: np.sqrt(k * (k + 2 * alpha - 1)
  1248. / (4 * (k + alpha) * (k + alpha - 1)))
  1249. f = lambda n, x: _ufuncs.eval_gegenbauer(n, alpha, x)
  1250. df = lambda n, x: ((-n*x*_ufuncs.eval_gegenbauer(n, alpha, x)
  1251. + ((n + 2*alpha - 1)
  1252. * _ufuncs.eval_gegenbauer(n - 1, alpha, x)))
  1253. / (1 - x**2))
  1254. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  1255. def gegenbauer(n, alpha, monic=False):
  1256. r"""Gegenbauer (ultraspherical) polynomial.
  1257. Defined to be the solution of
  1258. .. math::
  1259. (1 - x^2)\frac{d^2}{dx^2}C_n^{(\alpha)}
  1260. - (2\alpha + 1)x\frac{d}{dx}C_n^{(\alpha)}
  1261. + n(n + 2\alpha)C_n^{(\alpha)} = 0
  1262. for :math:`\alpha > -1/2`; :math:`C_n^{(\alpha)}` is a polynomial
  1263. of degree :math:`n`.
  1264. Parameters
  1265. ----------
  1266. n : int
  1267. Degree of the polynomial.
  1268. alpha : float
  1269. Parameter, must be greater than -0.5.
  1270. monic : bool, optional
  1271. If `True`, scale the leading coefficient to be 1. Default is
  1272. `False`.
  1273. Returns
  1274. -------
  1275. C : orthopoly1d
  1276. Gegenbauer polynomial.
  1277. Notes
  1278. -----
  1279. The polynomials :math:`C_n^{(\alpha)}` are orthogonal over
  1280. :math:`[-1,1]` with weight function :math:`(1 - x^2)^{(\alpha -
  1281. 1/2)}`.
  1282. Examples
  1283. --------
  1284. >>> import numpy as np
  1285. >>> from scipy import special
  1286. >>> import matplotlib.pyplot as plt
  1287. We can initialize a variable ``p`` as a Gegenbauer polynomial using the
  1288. `gegenbauer` function and evaluate at a point ``x = 1``.
  1289. >>> p = special.gegenbauer(3, 0.5, monic=False)
  1290. >>> p
  1291. poly1d([ 2.5, 0. , -1.5, 0. ])
  1292. >>> p(1)
  1293. 1.0
  1294. To evaluate ``p`` at various points ``x`` in the interval ``(-3, 3)``,
  1295. simply pass an array ``x`` to ``p`` as follows:
  1296. >>> x = np.linspace(-3, 3, 400)
  1297. >>> y = p(x)
  1298. We can then visualize ``x, y`` using `matplotlib.pyplot`.
  1299. >>> fig, ax = plt.subplots()
  1300. >>> ax.plot(x, y)
  1301. >>> ax.set_title("Gegenbauer (ultraspherical) polynomial of degree 3")
  1302. >>> ax.set_xlabel("x")
  1303. >>> ax.set_ylabel("G_3(x)")
  1304. >>> plt.show()
  1305. """
  1306. base = jacobi(n, alpha - 0.5, alpha - 0.5, monic=monic)
  1307. if monic:
  1308. return base
  1309. # Abrahmowitz and Stegan 22.5.20
  1310. factor = (_gam(2*alpha + n) * _gam(alpha + 0.5) /
  1311. _gam(2*alpha) / _gam(alpha + 0.5 + n))
  1312. base._scale(factor)
  1313. base.__dict__['_eval_func'] = lambda x: _ufuncs.eval_gegenbauer(float(n),
  1314. alpha, x)
  1315. return base
  1316. # Chebyshev of the first kind: T_n(x) =
  1317. # n! sqrt(pi) / _gam(n+1./2)* P^(-1/2,-1/2)_n(x)
  1318. # Computed anew.
  1319. def roots_chebyt(n, mu=False):
  1320. r"""Gauss-Chebyshev (first kind) quadrature.
  1321. Computes the sample points and weights for Gauss-Chebyshev
  1322. quadrature. The sample points are the roots of the nth degree
  1323. Chebyshev polynomial of the first kind, :math:`T_n(x)`. These
  1324. sample points and weights correctly integrate polynomials of
  1325. degree :math:`2n - 1` or less over the interval :math:`[-1, 1]`
  1326. with weight function :math:`w(x) = 1/\sqrt{1 - x^2}`. See 22.2.4
  1327. in [AS]_ for more details.
  1328. Parameters
  1329. ----------
  1330. n : int
  1331. quadrature order
  1332. mu : bool, optional
  1333. If True, return the sum of the weights, optional.
  1334. Returns
  1335. -------
  1336. x : ndarray
  1337. Sample points
  1338. w : ndarray
  1339. Weights
  1340. mu : float
  1341. Sum of the weights
  1342. See Also
  1343. --------
  1344. scipy.integrate.quadrature
  1345. scipy.integrate.fixed_quad
  1346. numpy.polynomial.chebyshev.chebgauss
  1347. References
  1348. ----------
  1349. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1350. Handbook of Mathematical Functions with Formulas,
  1351. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1352. """
  1353. m = int(n)
  1354. if n < 1 or n != m:
  1355. raise ValueError('n must be a positive integer.')
  1356. x = _ufuncs._sinpi(np.arange(-m + 1, m, 2) / (2*m))
  1357. w = np.full_like(x, pi/m)
  1358. if mu:
  1359. return x, w, pi
  1360. else:
  1361. return x, w
  1362. def chebyt(n, monic=False):
  1363. r"""Chebyshev polynomial of the first kind.
  1364. Defined to be the solution of
  1365. .. math::
  1366. (1 - x^2)\frac{d^2}{dx^2}T_n - x\frac{d}{dx}T_n + n^2T_n = 0;
  1367. :math:`T_n` is a polynomial of degree :math:`n`.
  1368. Parameters
  1369. ----------
  1370. n : int
  1371. Degree of the polynomial.
  1372. monic : bool, optional
  1373. If `True`, scale the leading coefficient to be 1. Default is
  1374. `False`.
  1375. Returns
  1376. -------
  1377. T : orthopoly1d
  1378. Chebyshev polynomial of the first kind.
  1379. Notes
  1380. -----
  1381. The polynomials :math:`T_n` are orthogonal over :math:`[-1, 1]`
  1382. with weight function :math:`(1 - x^2)^{-1/2}`.
  1383. See Also
  1384. --------
  1385. chebyu : Chebyshev polynomial of the second kind.
  1386. References
  1387. ----------
  1388. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1389. Handbook of Mathematical Functions with Formulas,
  1390. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1391. Examples
  1392. --------
  1393. Chebyshev polynomials of the first kind of order :math:`n` can
  1394. be obtained as the determinant of specific :math:`n \times n`
  1395. matrices. As an example we can check how the points obtained from
  1396. the determinant of the following :math:`3 \times 3` matrix
  1397. lay exacty on :math:`T_3`:
  1398. >>> import numpy as np
  1399. >>> import matplotlib.pyplot as plt
  1400. >>> from scipy.linalg import det
  1401. >>> from scipy.special import chebyt
  1402. >>> x = np.arange(-1.0, 1.0, 0.01)
  1403. >>> fig, ax = plt.subplots()
  1404. >>> ax.set_ylim(-2.0, 2.0)
  1405. >>> ax.set_title(r'Chebyshev polynomial $T_3$')
  1406. >>> ax.plot(x, chebyt(3)(x), label=rf'$T_3$')
  1407. >>> for p in np.arange(-1.0, 1.0, 0.1):
  1408. ... ax.plot(p,
  1409. ... det(np.array([[p, 1, 0], [1, 2*p, 1], [0, 1, 2*p]])),
  1410. ... 'rx')
  1411. >>> plt.legend(loc='best')
  1412. >>> plt.show()
  1413. They are also related to the Jacobi Polynomials
  1414. :math:`P_n^{(-0.5, -0.5)}` through the relation:
  1415. .. math::
  1416. P_n^{(-0.5, -0.5)}(x) = \frac{1}{4^n} \binom{2n}{n} T_n(x)
  1417. Let's verify it for :math:`n = 3`:
  1418. >>> from scipy.special import binom
  1419. >>> from scipy.special import jacobi
  1420. >>> x = np.arange(-1.0, 1.0, 0.01)
  1421. >>> np.allclose(jacobi(3, -0.5, -0.5)(x),
  1422. ... 1/64 * binom(6, 3) * chebyt(3)(x))
  1423. True
  1424. We can plot the Chebyshev polynomials :math:`T_n` for some values
  1425. of :math:`n`:
  1426. >>> x = np.arange(-1.5, 1.5, 0.01)
  1427. >>> fig, ax = plt.subplots()
  1428. >>> ax.set_ylim(-4.0, 4.0)
  1429. >>> ax.set_title(r'Chebyshev polynomials $T_n$')
  1430. >>> for n in np.arange(2,5):
  1431. ... ax.plot(x, chebyt(n)(x), label=rf'$T_n={n}$')
  1432. >>> plt.legend(loc='best')
  1433. >>> plt.show()
  1434. """
  1435. if n < 0:
  1436. raise ValueError("n must be nonnegative.")
  1437. wfunc = lambda x: 1.0 / sqrt(1 - x * x)
  1438. if n == 0:
  1439. return orthopoly1d([], [], pi, 1.0, wfunc, (-1, 1), monic,
  1440. lambda x: _ufuncs.eval_chebyt(n, x))
  1441. n1 = n
  1442. x, w, mu = roots_chebyt(n1, mu=True)
  1443. hn = pi / 2
  1444. kn = 2**(n - 1)
  1445. p = orthopoly1d(x, w, hn, kn, wfunc, (-1, 1), monic,
  1446. lambda x: _ufuncs.eval_chebyt(n, x))
  1447. return p
  1448. # Chebyshev of the second kind
  1449. # U_n(x) = (n+1)! sqrt(pi) / (2*_gam(n+3./2)) * P^(1/2,1/2)_n(x)
  1450. def roots_chebyu(n, mu=False):
  1451. r"""Gauss-Chebyshev (second kind) quadrature.
  1452. Computes the sample points and weights for Gauss-Chebyshev
  1453. quadrature. The sample points are the roots of the nth degree
  1454. Chebyshev polynomial of the second kind, :math:`U_n(x)`. These
  1455. sample points and weights correctly integrate polynomials of
  1456. degree :math:`2n - 1` or less over the interval :math:`[-1, 1]`
  1457. with weight function :math:`w(x) = \sqrt{1 - x^2}`. See 22.2.5 in
  1458. [AS]_ for details.
  1459. Parameters
  1460. ----------
  1461. n : int
  1462. quadrature order
  1463. mu : bool, optional
  1464. If True, return the sum of the weights, optional.
  1465. Returns
  1466. -------
  1467. x : ndarray
  1468. Sample points
  1469. w : ndarray
  1470. Weights
  1471. mu : float
  1472. Sum of the weights
  1473. See Also
  1474. --------
  1475. scipy.integrate.quadrature
  1476. scipy.integrate.fixed_quad
  1477. References
  1478. ----------
  1479. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1480. Handbook of Mathematical Functions with Formulas,
  1481. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1482. """
  1483. m = int(n)
  1484. if n < 1 or n != m:
  1485. raise ValueError('n must be a positive integer.')
  1486. t = np.arange(m, 0, -1) * pi / (m + 1)
  1487. x = np.cos(t)
  1488. w = pi * np.sin(t)**2 / (m + 1)
  1489. if mu:
  1490. return x, w, pi / 2
  1491. else:
  1492. return x, w
  1493. def chebyu(n, monic=False):
  1494. r"""Chebyshev polynomial of the second kind.
  1495. Defined to be the solution of
  1496. .. math::
  1497. (1 - x^2)\frac{d^2}{dx^2}U_n - 3x\frac{d}{dx}U_n
  1498. + n(n + 2)U_n = 0;
  1499. :math:`U_n` is a polynomial of degree :math:`n`.
  1500. Parameters
  1501. ----------
  1502. n : int
  1503. Degree of the polynomial.
  1504. monic : bool, optional
  1505. If `True`, scale the leading coefficient to be 1. Default is
  1506. `False`.
  1507. Returns
  1508. -------
  1509. U : orthopoly1d
  1510. Chebyshev polynomial of the second kind.
  1511. Notes
  1512. -----
  1513. The polynomials :math:`U_n` are orthogonal over :math:`[-1, 1]`
  1514. with weight function :math:`(1 - x^2)^{1/2}`.
  1515. See Also
  1516. --------
  1517. chebyt : Chebyshev polynomial of the first kind.
  1518. References
  1519. ----------
  1520. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1521. Handbook of Mathematical Functions with Formulas,
  1522. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1523. Examples
  1524. --------
  1525. Chebyshev polynomials of the second kind of order :math:`n` can
  1526. be obtained as the determinant of specific :math:`n \times n`
  1527. matrices. As an example we can check how the points obtained from
  1528. the determinant of the following :math:`3 \times 3` matrix
  1529. lay exacty on :math:`U_3`:
  1530. >>> import numpy as np
  1531. >>> import matplotlib.pyplot as plt
  1532. >>> from scipy.linalg import det
  1533. >>> from scipy.special import chebyu
  1534. >>> x = np.arange(-1.0, 1.0, 0.01)
  1535. >>> fig, ax = plt.subplots()
  1536. >>> ax.set_ylim(-2.0, 2.0)
  1537. >>> ax.set_title(r'Chebyshev polynomial $U_3$')
  1538. >>> ax.plot(x, chebyu(3)(x), label=rf'$U_3$')
  1539. >>> for p in np.arange(-1.0, 1.0, 0.1):
  1540. ... ax.plot(p,
  1541. ... det(np.array([[2*p, 1, 0], [1, 2*p, 1], [0, 1, 2*p]])),
  1542. ... 'rx')
  1543. >>> plt.legend(loc='best')
  1544. >>> plt.show()
  1545. They satisfy the recurrence relation:
  1546. .. math::
  1547. U_{2n-1}(x) = 2 T_n(x)U_{n-1}(x)
  1548. where the :math:`T_n` are the Chebyshev polynomial of the first kind.
  1549. Let's verify it for :math:`n = 2`:
  1550. >>> from scipy.special import chebyt
  1551. >>> x = np.arange(-1.0, 1.0, 0.01)
  1552. >>> np.allclose(chebyu(3)(x), 2 * chebyt(2)(x) * chebyu(1)(x))
  1553. True
  1554. We can plot the Chebyshev polynomials :math:`U_n` for some values
  1555. of :math:`n`:
  1556. >>> x = np.arange(-1.0, 1.0, 0.01)
  1557. >>> fig, ax = plt.subplots()
  1558. >>> ax.set_ylim(-1.5, 1.5)
  1559. >>> ax.set_title(r'Chebyshev polynomials $U_n$')
  1560. >>> for n in np.arange(1,5):
  1561. ... ax.plot(x, chebyu(n)(x), label=rf'$U_n={n}$')
  1562. >>> plt.legend(loc='best')
  1563. >>> plt.show()
  1564. """
  1565. base = jacobi(n, 0.5, 0.5, monic=monic)
  1566. if monic:
  1567. return base
  1568. factor = sqrt(pi) / 2.0 * _gam(n + 2) / _gam(n + 1.5)
  1569. base._scale(factor)
  1570. return base
  1571. # Chebyshev of the first kind C_n(x)
  1572. def roots_chebyc(n, mu=False):
  1573. r"""Gauss-Chebyshev (first kind) quadrature.
  1574. Compute the sample points and weights for Gauss-Chebyshev
  1575. quadrature. The sample points are the roots of the nth degree
  1576. Chebyshev polynomial of the first kind, :math:`C_n(x)`. These
  1577. sample points and weights correctly integrate polynomials of
  1578. degree :math:`2n - 1` or less over the interval :math:`[-2, 2]`
  1579. with weight function :math:`w(x) = 1 / \sqrt{1 - (x/2)^2}`. See
  1580. 22.2.6 in [AS]_ for more details.
  1581. Parameters
  1582. ----------
  1583. n : int
  1584. quadrature order
  1585. mu : bool, optional
  1586. If True, return the sum of the weights, optional.
  1587. Returns
  1588. -------
  1589. x : ndarray
  1590. Sample points
  1591. w : ndarray
  1592. Weights
  1593. mu : float
  1594. Sum of the weights
  1595. See Also
  1596. --------
  1597. scipy.integrate.quadrature
  1598. scipy.integrate.fixed_quad
  1599. References
  1600. ----------
  1601. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1602. Handbook of Mathematical Functions with Formulas,
  1603. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1604. """
  1605. x, w, m = roots_chebyt(n, True)
  1606. x *= 2
  1607. w *= 2
  1608. m *= 2
  1609. if mu:
  1610. return x, w, m
  1611. else:
  1612. return x, w
  1613. def chebyc(n, monic=False):
  1614. r"""Chebyshev polynomial of the first kind on :math:`[-2, 2]`.
  1615. Defined as :math:`C_n(x) = 2T_n(x/2)`, where :math:`T_n` is the
  1616. nth Chebychev polynomial of the first kind.
  1617. Parameters
  1618. ----------
  1619. n : int
  1620. Degree of the polynomial.
  1621. monic : bool, optional
  1622. If `True`, scale the leading coefficient to be 1. Default is
  1623. `False`.
  1624. Returns
  1625. -------
  1626. C : orthopoly1d
  1627. Chebyshev polynomial of the first kind on :math:`[-2, 2]`.
  1628. Notes
  1629. -----
  1630. The polynomials :math:`C_n(x)` are orthogonal over :math:`[-2, 2]`
  1631. with weight function :math:`1/\sqrt{1 - (x/2)^2}`.
  1632. See Also
  1633. --------
  1634. chebyt : Chebyshev polynomial of the first kind.
  1635. References
  1636. ----------
  1637. .. [1] Abramowitz and Stegun, "Handbook of Mathematical Functions"
  1638. Section 22. National Bureau of Standards, 1972.
  1639. """
  1640. if n < 0:
  1641. raise ValueError("n must be nonnegative.")
  1642. if n == 0:
  1643. n1 = n + 1
  1644. else:
  1645. n1 = n
  1646. x, w = roots_chebyc(n1)
  1647. if n == 0:
  1648. x, w = [], []
  1649. hn = 4 * pi * ((n == 0) + 1)
  1650. kn = 1.0
  1651. p = orthopoly1d(x, w, hn, kn,
  1652. wfunc=lambda x: 1.0 / sqrt(1 - x * x / 4.0),
  1653. limits=(-2, 2), monic=monic)
  1654. if not monic:
  1655. p._scale(2.0 / p(2))
  1656. p.__dict__['_eval_func'] = lambda x: _ufuncs.eval_chebyc(n, x)
  1657. return p
  1658. # Chebyshev of the second kind S_n(x)
  1659. def roots_chebys(n, mu=False):
  1660. r"""Gauss-Chebyshev (second kind) quadrature.
  1661. Compute the sample points and weights for Gauss-Chebyshev
  1662. quadrature. The sample points are the roots of the nth degree
  1663. Chebyshev polynomial of the second kind, :math:`S_n(x)`. These
  1664. sample points and weights correctly integrate polynomials of
  1665. degree :math:`2n - 1` or less over the interval :math:`[-2, 2]`
  1666. with weight function :math:`w(x) = \sqrt{1 - (x/2)^2}`. See 22.2.7
  1667. in [AS]_ for more details.
  1668. Parameters
  1669. ----------
  1670. n : int
  1671. quadrature order
  1672. mu : bool, optional
  1673. If True, return the sum of the weights, optional.
  1674. Returns
  1675. -------
  1676. x : ndarray
  1677. Sample points
  1678. w : ndarray
  1679. Weights
  1680. mu : float
  1681. Sum of the weights
  1682. See Also
  1683. --------
  1684. scipy.integrate.quadrature
  1685. scipy.integrate.fixed_quad
  1686. References
  1687. ----------
  1688. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1689. Handbook of Mathematical Functions with Formulas,
  1690. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1691. """
  1692. x, w, m = roots_chebyu(n, True)
  1693. x *= 2
  1694. w *= 2
  1695. m *= 2
  1696. if mu:
  1697. return x, w, m
  1698. else:
  1699. return x, w
  1700. def chebys(n, monic=False):
  1701. r"""Chebyshev polynomial of the second kind on :math:`[-2, 2]`.
  1702. Defined as :math:`S_n(x) = U_n(x/2)` where :math:`U_n` is the
  1703. nth Chebychev polynomial of the second kind.
  1704. Parameters
  1705. ----------
  1706. n : int
  1707. Degree of the polynomial.
  1708. monic : bool, optional
  1709. If `True`, scale the leading coefficient to be 1. Default is
  1710. `False`.
  1711. Returns
  1712. -------
  1713. S : orthopoly1d
  1714. Chebyshev polynomial of the second kind on :math:`[-2, 2]`.
  1715. Notes
  1716. -----
  1717. The polynomials :math:`S_n(x)` are orthogonal over :math:`[-2, 2]`
  1718. with weight function :math:`\sqrt{1 - (x/2)}^2`.
  1719. See Also
  1720. --------
  1721. chebyu : Chebyshev polynomial of the second kind
  1722. References
  1723. ----------
  1724. .. [1] Abramowitz and Stegun, "Handbook of Mathematical Functions"
  1725. Section 22. National Bureau of Standards, 1972.
  1726. """
  1727. if n < 0:
  1728. raise ValueError("n must be nonnegative.")
  1729. if n == 0:
  1730. n1 = n + 1
  1731. else:
  1732. n1 = n
  1733. x, w = roots_chebys(n1)
  1734. if n == 0:
  1735. x, w = [], []
  1736. hn = pi
  1737. kn = 1.0
  1738. p = orthopoly1d(x, w, hn, kn,
  1739. wfunc=lambda x: sqrt(1 - x * x / 4.0),
  1740. limits=(-2, 2), monic=monic)
  1741. if not monic:
  1742. factor = (n + 1.0) / p(2)
  1743. p._scale(factor)
  1744. p.__dict__['_eval_func'] = lambda x: _ufuncs.eval_chebys(n, x)
  1745. return p
  1746. # Shifted Chebyshev of the first kind T^*_n(x)
  1747. def roots_sh_chebyt(n, mu=False):
  1748. r"""Gauss-Chebyshev (first kind, shifted) quadrature.
  1749. Compute the sample points and weights for Gauss-Chebyshev
  1750. quadrature. The sample points are the roots of the nth degree
  1751. shifted Chebyshev polynomial of the first kind, :math:`T_n(x)`.
  1752. These sample points and weights correctly integrate polynomials of
  1753. degree :math:`2n - 1` or less over the interval :math:`[0, 1]`
  1754. with weight function :math:`w(x) = 1/\sqrt{x - x^2}`. See 22.2.8
  1755. in [AS]_ for more details.
  1756. Parameters
  1757. ----------
  1758. n : int
  1759. quadrature order
  1760. mu : bool, optional
  1761. If True, return the sum of the weights, optional.
  1762. Returns
  1763. -------
  1764. x : ndarray
  1765. Sample points
  1766. w : ndarray
  1767. Weights
  1768. mu : float
  1769. Sum of the weights
  1770. See Also
  1771. --------
  1772. scipy.integrate.quadrature
  1773. scipy.integrate.fixed_quad
  1774. References
  1775. ----------
  1776. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1777. Handbook of Mathematical Functions with Formulas,
  1778. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1779. """
  1780. xw = roots_chebyt(n, mu)
  1781. return ((xw[0] + 1) / 2,) + xw[1:]
  1782. def sh_chebyt(n, monic=False):
  1783. r"""Shifted Chebyshev polynomial of the first kind.
  1784. Defined as :math:`T^*_n(x) = T_n(2x - 1)` for :math:`T_n` the nth
  1785. Chebyshev polynomial of the first kind.
  1786. Parameters
  1787. ----------
  1788. n : int
  1789. Degree of the polynomial.
  1790. monic : bool, optional
  1791. If `True`, scale the leading coefficient to be 1. Default is
  1792. `False`.
  1793. Returns
  1794. -------
  1795. T : orthopoly1d
  1796. Shifted Chebyshev polynomial of the first kind.
  1797. Notes
  1798. -----
  1799. The polynomials :math:`T^*_n` are orthogonal over :math:`[0, 1]`
  1800. with weight function :math:`(x - x^2)^{-1/2}`.
  1801. """
  1802. base = sh_jacobi(n, 0.0, 0.5, monic=monic)
  1803. if monic:
  1804. return base
  1805. if n > 0:
  1806. factor = 4**n / 2.0
  1807. else:
  1808. factor = 1.0
  1809. base._scale(factor)
  1810. return base
  1811. # Shifted Chebyshev of the second kind U^*_n(x)
  1812. def roots_sh_chebyu(n, mu=False):
  1813. r"""Gauss-Chebyshev (second kind, shifted) quadrature.
  1814. Computes the sample points and weights for Gauss-Chebyshev
  1815. quadrature. The sample points are the roots of the nth degree
  1816. shifted Chebyshev polynomial of the second kind, :math:`U_n(x)`.
  1817. These sample points and weights correctly integrate polynomials of
  1818. degree :math:`2n - 1` or less over the interval :math:`[0, 1]`
  1819. with weight function :math:`w(x) = \sqrt{x - x^2}`. See 22.2.9 in
  1820. [AS]_ for more details.
  1821. Parameters
  1822. ----------
  1823. n : int
  1824. quadrature order
  1825. mu : bool, optional
  1826. If True, return the sum of the weights, optional.
  1827. Returns
  1828. -------
  1829. x : ndarray
  1830. Sample points
  1831. w : ndarray
  1832. Weights
  1833. mu : float
  1834. Sum of the weights
  1835. See Also
  1836. --------
  1837. scipy.integrate.quadrature
  1838. scipy.integrate.fixed_quad
  1839. References
  1840. ----------
  1841. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1842. Handbook of Mathematical Functions with Formulas,
  1843. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1844. """
  1845. x, w, m = roots_chebyu(n, True)
  1846. x = (x + 1) / 2
  1847. m_us = _ufuncs.beta(1.5, 1.5)
  1848. w *= m_us / m
  1849. if mu:
  1850. return x, w, m_us
  1851. else:
  1852. return x, w
  1853. def sh_chebyu(n, monic=False):
  1854. r"""Shifted Chebyshev polynomial of the second kind.
  1855. Defined as :math:`U^*_n(x) = U_n(2x - 1)` for :math:`U_n` the nth
  1856. Chebyshev polynomial of the second kind.
  1857. Parameters
  1858. ----------
  1859. n : int
  1860. Degree of the polynomial.
  1861. monic : bool, optional
  1862. If `True`, scale the leading coefficient to be 1. Default is
  1863. `False`.
  1864. Returns
  1865. -------
  1866. U : orthopoly1d
  1867. Shifted Chebyshev polynomial of the second kind.
  1868. Notes
  1869. -----
  1870. The polynomials :math:`U^*_n` are orthogonal over :math:`[0, 1]`
  1871. with weight function :math:`(x - x^2)^{1/2}`.
  1872. """
  1873. base = sh_jacobi(n, 2.0, 1.5, monic=monic)
  1874. if monic:
  1875. return base
  1876. factor = 4**n
  1877. base._scale(factor)
  1878. return base
  1879. # Legendre
  1880. def roots_legendre(n, mu=False):
  1881. r"""Gauss-Legendre quadrature.
  1882. Compute the sample points and weights for Gauss-Legendre
  1883. quadrature [GL]_. The sample points are the roots of the nth degree
  1884. Legendre polynomial :math:`P_n(x)`. These sample points and
  1885. weights correctly integrate polynomials of degree :math:`2n - 1`
  1886. or less over the interval :math:`[-1, 1]` with weight function
  1887. :math:`w(x) = 1`. See 2.2.10 in [AS]_ for more details.
  1888. Parameters
  1889. ----------
  1890. n : int
  1891. quadrature order
  1892. mu : bool, optional
  1893. If True, return the sum of the weights, optional.
  1894. Returns
  1895. -------
  1896. x : ndarray
  1897. Sample points
  1898. w : ndarray
  1899. Weights
  1900. mu : float
  1901. Sum of the weights
  1902. See Also
  1903. --------
  1904. scipy.integrate.quadrature
  1905. scipy.integrate.fixed_quad
  1906. numpy.polynomial.legendre.leggauss
  1907. References
  1908. ----------
  1909. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  1910. Handbook of Mathematical Functions with Formulas,
  1911. Graphs, and Mathematical Tables. New York: Dover, 1972.
  1912. .. [GL] Gauss-Legendre quadrature, Wikipedia,
  1913. https://en.wikipedia.org/wiki/Gauss%E2%80%93Legendre_quadrature
  1914. Examples
  1915. --------
  1916. >>> import numpy as np
  1917. >>> from scipy.special import roots_legendre, eval_legendre
  1918. >>> roots, weights = roots_legendre(9)
  1919. ``roots`` holds the roots, and ``weights`` holds the weights for
  1920. Gauss-Legendre quadrature.
  1921. >>> roots
  1922. array([-0.96816024, -0.83603111, -0.61337143, -0.32425342, 0. ,
  1923. 0.32425342, 0.61337143, 0.83603111, 0.96816024])
  1924. >>> weights
  1925. array([0.08127439, 0.18064816, 0.2606107 , 0.31234708, 0.33023936,
  1926. 0.31234708, 0.2606107 , 0.18064816, 0.08127439])
  1927. Verify that we have the roots by evaluating the degree 9 Legendre
  1928. polynomial at ``roots``. All the values are approximately zero:
  1929. >>> eval_legendre(9, roots)
  1930. array([-8.88178420e-16, -2.22044605e-16, 1.11022302e-16, 1.11022302e-16,
  1931. 0.00000000e+00, -5.55111512e-17, -1.94289029e-16, 1.38777878e-16,
  1932. -8.32667268e-17])
  1933. Here we'll show how the above values can be used to estimate the
  1934. integral from 1 to 2 of f(t) = t + 1/t with Gauss-Legendre
  1935. quadrature [GL]_. First define the function and the integration
  1936. limits.
  1937. >>> def f(t):
  1938. ... return t + 1/t
  1939. ...
  1940. >>> a = 1
  1941. >>> b = 2
  1942. We'll use ``integral(f(t), t=a, t=b)`` to denote the definite integral
  1943. of f from t=a to t=b. The sample points in ``roots`` are from the
  1944. interval [-1, 1], so we'll rewrite the integral with the simple change
  1945. of variable::
  1946. x = 2/(b - a) * t - (a + b)/(b - a)
  1947. with inverse::
  1948. t = (b - a)/2 * x + (a + 2)/2
  1949. Then::
  1950. integral(f(t), a, b) =
  1951. (b - a)/2 * integral(f((b-a)/2*x + (a+b)/2), x=-1, x=1)
  1952. We can approximate the latter integral with the values returned
  1953. by `roots_legendre`.
  1954. Map the roots computed above from [-1, 1] to [a, b].
  1955. >>> t = (b - a)/2 * roots + (a + b)/2
  1956. Approximate the integral as the weighted sum of the function values.
  1957. >>> (b - a)/2 * f(t).dot(weights)
  1958. 2.1931471805599276
  1959. Compare that to the exact result, which is 3/2 + log(2):
  1960. >>> 1.5 + np.log(2)
  1961. 2.1931471805599454
  1962. """
  1963. m = int(n)
  1964. if n < 1 or n != m:
  1965. raise ValueError("n must be a positive integer.")
  1966. mu0 = 2.0
  1967. an_func = lambda k: 0.0 * k
  1968. bn_func = lambda k: k * np.sqrt(1.0 / (4 * k * k - 1))
  1969. f = _ufuncs.eval_legendre
  1970. df = lambda n, x: (-n*x*_ufuncs.eval_legendre(n, x)
  1971. + n*_ufuncs.eval_legendre(n-1, x))/(1-x**2)
  1972. return _gen_roots_and_weights(m, mu0, an_func, bn_func, f, df, True, mu)
  1973. def legendre(n, monic=False):
  1974. r"""Legendre polynomial.
  1975. Defined to be the solution of
  1976. .. math::
  1977. \frac{d}{dx}\left[(1 - x^2)\frac{d}{dx}P_n(x)\right]
  1978. + n(n + 1)P_n(x) = 0;
  1979. :math:`P_n(x)` is a polynomial of degree :math:`n`.
  1980. Parameters
  1981. ----------
  1982. n : int
  1983. Degree of the polynomial.
  1984. monic : bool, optional
  1985. If `True`, scale the leading coefficient to be 1. Default is
  1986. `False`.
  1987. Returns
  1988. -------
  1989. P : orthopoly1d
  1990. Legendre polynomial.
  1991. Notes
  1992. -----
  1993. The polynomials :math:`P_n` are orthogonal over :math:`[-1, 1]`
  1994. with weight function 1.
  1995. Examples
  1996. --------
  1997. Generate the 3rd-order Legendre polynomial 1/2*(5x^3 + 0x^2 - 3x + 0):
  1998. >>> from scipy.special import legendre
  1999. >>> legendre(3)
  2000. poly1d([ 2.5, 0. , -1.5, 0. ])
  2001. """
  2002. if n < 0:
  2003. raise ValueError("n must be nonnegative.")
  2004. if n == 0:
  2005. n1 = n + 1
  2006. else:
  2007. n1 = n
  2008. x, w = roots_legendre(n1)
  2009. if n == 0:
  2010. x, w = [], []
  2011. hn = 2.0 / (2 * n + 1)
  2012. kn = _gam(2 * n + 1) / _gam(n + 1)**2 / 2.0**n
  2013. p = orthopoly1d(x, w, hn, kn, wfunc=lambda x: 1.0, limits=(-1, 1),
  2014. monic=monic,
  2015. eval_func=lambda x: _ufuncs.eval_legendre(n, x))
  2016. return p
  2017. # Shifted Legendre P^*_n(x)
  2018. def roots_sh_legendre(n, mu=False):
  2019. r"""Gauss-Legendre (shifted) quadrature.
  2020. Compute the sample points and weights for Gauss-Legendre
  2021. quadrature. The sample points are the roots of the nth degree
  2022. shifted Legendre polynomial :math:`P^*_n(x)`. These sample points
  2023. and weights correctly integrate polynomials of degree :math:`2n -
  2024. 1` or less over the interval :math:`[0, 1]` with weight function
  2025. :math:`w(x) = 1.0`. See 2.2.11 in [AS]_ for details.
  2026. Parameters
  2027. ----------
  2028. n : int
  2029. quadrature order
  2030. mu : bool, optional
  2031. If True, return the sum of the weights, optional.
  2032. Returns
  2033. -------
  2034. x : ndarray
  2035. Sample points
  2036. w : ndarray
  2037. Weights
  2038. mu : float
  2039. Sum of the weights
  2040. See Also
  2041. --------
  2042. scipy.integrate.quadrature
  2043. scipy.integrate.fixed_quad
  2044. References
  2045. ----------
  2046. .. [AS] Milton Abramowitz and Irene A. Stegun, eds.
  2047. Handbook of Mathematical Functions with Formulas,
  2048. Graphs, and Mathematical Tables. New York: Dover, 1972.
  2049. """
  2050. x, w = roots_legendre(n)
  2051. x = (x + 1) / 2
  2052. w /= 2
  2053. if mu:
  2054. return x, w, 1.0
  2055. else:
  2056. return x, w
  2057. def sh_legendre(n, monic=False):
  2058. r"""Shifted Legendre polynomial.
  2059. Defined as :math:`P^*_n(x) = P_n(2x - 1)` for :math:`P_n` the nth
  2060. Legendre polynomial.
  2061. Parameters
  2062. ----------
  2063. n : int
  2064. Degree of the polynomial.
  2065. monic : bool, optional
  2066. If `True`, scale the leading coefficient to be 1. Default is
  2067. `False`.
  2068. Returns
  2069. -------
  2070. P : orthopoly1d
  2071. Shifted Legendre polynomial.
  2072. Notes
  2073. -----
  2074. The polynomials :math:`P^*_n` are orthogonal over :math:`[0, 1]`
  2075. with weight function 1.
  2076. """
  2077. if n < 0:
  2078. raise ValueError("n must be nonnegative.")
  2079. wfunc = lambda x: 0.0 * x + 1.0
  2080. if n == 0:
  2081. return orthopoly1d([], [], 1.0, 1.0, wfunc, (0, 1), monic,
  2082. lambda x: _ufuncs.eval_sh_legendre(n, x))
  2083. x, w = roots_sh_legendre(n)
  2084. hn = 1.0 / (2 * n + 1.0)
  2085. kn = _gam(2 * n + 1) / _gam(n + 1)**2
  2086. p = orthopoly1d(x, w, hn, kn, wfunc, limits=(0, 1), monic=monic,
  2087. eval_func=lambda x: _ufuncs.eval_sh_legendre(n, x))
  2088. return p
  2089. # Make the old root function names an alias for the new ones
  2090. _modattrs = globals()
  2091. for newfun, oldfun in _rootfuns_map.items():
  2092. _modattrs[oldfun] = _modattrs[newfun]
  2093. __all__.append(oldfun)