chebyshev.py 61 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807180818091810181118121813181418151816181718181819182018211822182318241825182618271828182918301831183218331834183518361837183818391840184118421843184418451846184718481849185018511852185318541855185618571858185918601861186218631864186518661867186818691870187118721873187418751876187718781879188018811882188318841885188618871888188918901891189218931894189518961897189818991900190119021903190419051906190719081909191019111912191319141915191619171918191919201921192219231924192519261927192819291930193119321933193419351936193719381939194019411942194319441945194619471948194919501951195219531954195519561957195819591960196119621963196419651966196719681969197019711972197319741975197619771978197919801981198219831984198519861987198819891990199119921993199419951996199719981999200020012002200320042005200620072008200920102011201220132014201520162017201820192020202120222023202420252026202720282029203020312032203320342035203620372038203920402041204220432044204520462047204820492050205120522053205420552056205720582059206020612062206320642065206620672068206920702071207220732074207520762077207820792080
  1. """
  2. Objects for dealing with Chebyshev series.
  3. This module provides a number of objects (mostly functions) useful for
  4. dealing with Chebyshev series, including a `Chebyshev` class that
  5. encapsulates the usual arithmetic operations. (General information
  6. on how this module represents and works with such polynomials is in the
  7. docstring for its "parent" sub-package, `numpy.polynomial`).
  8. Constants
  9. ---------
  10. - `chebdomain` -- Chebyshev series default domain, [-1,1].
  11. - `chebzero` -- (Coefficients of the) Chebyshev series that evaluates
  12. identically to 0.
  13. - `chebone` -- (Coefficients of the) Chebyshev series that evaluates
  14. identically to 1.
  15. - `chebx` -- (Coefficients of the) Chebyshev series for the identity map,
  16. ``f(x) = x``.
  17. Arithmetic
  18. ----------
  19. - `chebadd` -- add two Chebyshev series.
  20. - `chebsub` -- subtract one Chebyshev series from another.
  21. - `chebmul` -- multiply two Chebyshev series.
  22. - `chebdiv` -- divide one Chebyshev series by another.
  23. - `chebpow` -- raise a Chebyshev series to an positive integer power
  24. - `chebval` -- evaluate a Chebyshev series at given points.
  25. - `chebval2d` -- evaluate a 2D Chebyshev series at given points.
  26. - `chebval3d` -- evaluate a 3D Chebyshev series at given points.
  27. - `chebgrid2d` -- evaluate a 2D Chebyshev series on a Cartesian product.
  28. - `chebgrid3d` -- evaluate a 3D Chebyshev series on a Cartesian product.
  29. Calculus
  30. --------
  31. - `chebder` -- differentiate a Chebyshev series.
  32. - `chebint` -- integrate a Chebyshev series.
  33. Misc Functions
  34. --------------
  35. - `chebfromroots` -- create a Chebyshev series with specified roots.
  36. - `chebroots` -- find the roots of a Chebyshev series.
  37. - `chebvander` -- Vandermonde-like matrix for Chebyshev polynomials.
  38. - `chebvander2d` -- Vandermonde-like matrix for 2D power series.
  39. - `chebvander3d` -- Vandermonde-like matrix for 3D power series.
  40. - `chebgauss` -- Gauss-Chebyshev quadrature, points and weights.
  41. - `chebweight` -- Chebyshev weight function.
  42. - `chebcompanion` -- symmetrized companion matrix in Chebyshev form.
  43. - `chebfit` -- least-squares fit returning a Chebyshev series.
  44. - `chebpts1` -- Chebyshev points of the first kind.
  45. - `chebpts2` -- Chebyshev points of the second kind.
  46. - `chebtrim` -- trim leading coefficients from a Chebyshev series.
  47. - `chebline` -- Chebyshev series representing given straight line.
  48. - `cheb2poly` -- convert a Chebyshev series to a polynomial.
  49. - `poly2cheb` -- convert a polynomial to a Chebyshev series.
  50. Classes
  51. -------
  52. - `Chebyshev` -- A Chebyshev series class.
  53. See also
  54. --------
  55. `numpy.polynomial`
  56. Notes
  57. -----
  58. The implementations of multiplication, division, integration, and
  59. differentiation use the algebraic identities [1]_:
  60. .. math ::
  61. T_n(x) = \\frac{z^n + z^{-n}}{2} \\\\
  62. z\\frac{dx}{dz} = \\frac{z - z^{-1}}{2}.
  63. where
  64. .. math :: x = \\frac{z + z^{-1}}{2}.
  65. These identities allow a Chebyshev series to be expressed as a finite,
  66. symmetric Laurent series. In this module, this sort of Laurent series
  67. is referred to as a "z-series."
  68. References
  69. ----------
  70. .. [1] A. T. Benjamin, et al., "Combinatorial Trigonometry with Chebyshev
  71. Polynomials," *Journal of Statistical Planning and Inference 14*, 2008
  72. (preprint: http://www.math.hmc.edu/~benjamin/papers/CombTrig.pdf, pg. 4)
  73. """
  74. from __future__ import division, absolute_import, print_function
  75. import warnings
  76. import numpy as np
  77. import numpy.linalg as la
  78. from . import polyutils as pu
  79. from ._polybase import ABCPolyBase
  80. __all__ = [
  81. 'chebzero', 'chebone', 'chebx', 'chebdomain', 'chebline', 'chebadd',
  82. 'chebsub', 'chebmulx', 'chebmul', 'chebdiv', 'chebpow', 'chebval',
  83. 'chebder', 'chebint', 'cheb2poly', 'poly2cheb', 'chebfromroots',
  84. 'chebvander', 'chebfit', 'chebtrim', 'chebroots', 'chebpts1',
  85. 'chebpts2', 'Chebyshev', 'chebval2d', 'chebval3d', 'chebgrid2d',
  86. 'chebgrid3d', 'chebvander2d', 'chebvander3d', 'chebcompanion',
  87. 'chebgauss', 'chebweight']
  88. chebtrim = pu.trimcoef
  89. #
  90. # A collection of functions for manipulating z-series. These are private
  91. # functions and do minimal error checking.
  92. #
  93. def _cseries_to_zseries(c):
  94. """Covert Chebyshev series to z-series.
  95. Covert a Chebyshev series to the equivalent z-series. The result is
  96. never an empty array. The dtype of the return is the same as that of
  97. the input. No checks are run on the arguments as this routine is for
  98. internal use.
  99. Parameters
  100. ----------
  101. c : 1-D ndarray
  102. Chebyshev coefficients, ordered from low to high
  103. Returns
  104. -------
  105. zs : 1-D ndarray
  106. Odd length symmetric z-series, ordered from low to high.
  107. """
  108. n = c.size
  109. zs = np.zeros(2*n-1, dtype=c.dtype)
  110. zs[n-1:] = c/2
  111. return zs + zs[::-1]
  112. def _zseries_to_cseries(zs):
  113. """Covert z-series to a Chebyshev series.
  114. Covert a z series to the equivalent Chebyshev series. The result is
  115. never an empty array. The dtype of the return is the same as that of
  116. the input. No checks are run on the arguments as this routine is for
  117. internal use.
  118. Parameters
  119. ----------
  120. zs : 1-D ndarray
  121. Odd length symmetric z-series, ordered from low to high.
  122. Returns
  123. -------
  124. c : 1-D ndarray
  125. Chebyshev coefficients, ordered from low to high.
  126. """
  127. n = (zs.size + 1)//2
  128. c = zs[n-1:].copy()
  129. c[1:n] *= 2
  130. return c
  131. def _zseries_mul(z1, z2):
  132. """Multiply two z-series.
  133. Multiply two z-series to produce a z-series.
  134. Parameters
  135. ----------
  136. z1, z2 : 1-D ndarray
  137. The arrays must be 1-D but this is not checked.
  138. Returns
  139. -------
  140. product : 1-D ndarray
  141. The product z-series.
  142. Notes
  143. -----
  144. This is simply convolution. If symmetric/anti-symmetric z-series are
  145. denoted by S/A then the following rules apply:
  146. S*S, A*A -> S
  147. S*A, A*S -> A
  148. """
  149. return np.convolve(z1, z2)
  150. def _zseries_div(z1, z2):
  151. """Divide the first z-series by the second.
  152. Divide `z1` by `z2` and return the quotient and remainder as z-series.
  153. Warning: this implementation only applies when both z1 and z2 have the
  154. same symmetry, which is sufficient for present purposes.
  155. Parameters
  156. ----------
  157. z1, z2 : 1-D ndarray
  158. The arrays must be 1-D and have the same symmetry, but this is not
  159. checked.
  160. Returns
  161. -------
  162. (quotient, remainder) : 1-D ndarrays
  163. Quotient and remainder as z-series.
  164. Notes
  165. -----
  166. This is not the same as polynomial division on account of the desired form
  167. of the remainder. If symmetric/anti-symmetric z-series are denoted by S/A
  168. then the following rules apply:
  169. S/S -> S,S
  170. A/A -> S,A
  171. The restriction to types of the same symmetry could be fixed but seems like
  172. unneeded generality. There is no natural form for the remainder in the case
  173. where there is no symmetry.
  174. """
  175. z1 = z1.copy()
  176. z2 = z2.copy()
  177. len1 = len(z1)
  178. len2 = len(z2)
  179. if len2 == 1:
  180. z1 /= z2
  181. return z1, z1[:1]*0
  182. elif len1 < len2:
  183. return z1[:1]*0, z1
  184. else:
  185. dlen = len1 - len2
  186. scl = z2[0]
  187. z2 /= scl
  188. quo = np.empty(dlen + 1, dtype=z1.dtype)
  189. i = 0
  190. j = dlen
  191. while i < j:
  192. r = z1[i]
  193. quo[i] = z1[i]
  194. quo[dlen - i] = r
  195. tmp = r*z2
  196. z1[i:i+len2] -= tmp
  197. z1[j:j+len2] -= tmp
  198. i += 1
  199. j -= 1
  200. r = z1[i]
  201. quo[i] = r
  202. tmp = r*z2
  203. z1[i:i+len2] -= tmp
  204. quo /= scl
  205. rem = z1[i+1:i-1+len2].copy()
  206. return quo, rem
  207. def _zseries_der(zs):
  208. """Differentiate a z-series.
  209. The derivative is with respect to x, not z. This is achieved using the
  210. chain rule and the value of dx/dz given in the module notes.
  211. Parameters
  212. ----------
  213. zs : z-series
  214. The z-series to differentiate.
  215. Returns
  216. -------
  217. derivative : z-series
  218. The derivative
  219. Notes
  220. -----
  221. The zseries for x (ns) has been multiplied by two in order to avoid
  222. using floats that are incompatible with Decimal and likely other
  223. specialized scalar types. This scaling has been compensated by
  224. multiplying the value of zs by two also so that the two cancels in the
  225. division.
  226. """
  227. n = len(zs)//2
  228. ns = np.array([-1, 0, 1], dtype=zs.dtype)
  229. zs *= np.arange(-n, n+1)*2
  230. d, r = _zseries_div(zs, ns)
  231. return d
  232. def _zseries_int(zs):
  233. """Integrate a z-series.
  234. The integral is with respect to x, not z. This is achieved by a change
  235. of variable using dx/dz given in the module notes.
  236. Parameters
  237. ----------
  238. zs : z-series
  239. The z-series to integrate
  240. Returns
  241. -------
  242. integral : z-series
  243. The indefinite integral
  244. Notes
  245. -----
  246. The zseries for x (ns) has been multiplied by two in order to avoid
  247. using floats that are incompatible with Decimal and likely other
  248. specialized scalar types. This scaling has been compensated by
  249. dividing the resulting zs by two.
  250. """
  251. n = 1 + len(zs)//2
  252. ns = np.array([-1, 0, 1], dtype=zs.dtype)
  253. zs = _zseries_mul(zs, ns)
  254. div = np.arange(-n, n+1)*2
  255. zs[:n] /= div[:n]
  256. zs[n+1:] /= div[n+1:]
  257. zs[n] = 0
  258. return zs
  259. #
  260. # Chebyshev series functions
  261. #
  262. def poly2cheb(pol):
  263. """
  264. Convert a polynomial to a Chebyshev series.
  265. Convert an array representing the coefficients of a polynomial (relative
  266. to the "standard" basis) ordered from lowest degree to highest, to an
  267. array of the coefficients of the equivalent Chebyshev series, ordered
  268. from lowest to highest degree.
  269. Parameters
  270. ----------
  271. pol : array_like
  272. 1-D array containing the polynomial coefficients
  273. Returns
  274. -------
  275. c : ndarray
  276. 1-D array containing the coefficients of the equivalent Chebyshev
  277. series.
  278. See Also
  279. --------
  280. cheb2poly
  281. Notes
  282. -----
  283. The easy way to do conversions between polynomial basis sets
  284. is to use the convert method of a class instance.
  285. Examples
  286. --------
  287. >>> from numpy import polynomial as P
  288. >>> p = P.Polynomial(range(4))
  289. >>> p
  290. Polynomial([ 0., 1., 2., 3.], [-1., 1.])
  291. >>> c = p.convert(kind=P.Chebyshev)
  292. >>> c
  293. Chebyshev([ 1. , 3.25, 1. , 0.75], [-1., 1.])
  294. >>> P.poly2cheb(range(4))
  295. array([ 1. , 3.25, 1. , 0.75])
  296. """
  297. [pol] = pu.as_series([pol])
  298. deg = len(pol) - 1
  299. res = 0
  300. for i in range(deg, -1, -1):
  301. res = chebadd(chebmulx(res), pol[i])
  302. return res
  303. def cheb2poly(c):
  304. """
  305. Convert a Chebyshev series to a polynomial.
  306. Convert an array representing the coefficients of a Chebyshev series,
  307. ordered from lowest degree to highest, to an array of the coefficients
  308. of the equivalent polynomial (relative to the "standard" basis) ordered
  309. from lowest to highest degree.
  310. Parameters
  311. ----------
  312. c : array_like
  313. 1-D array containing the Chebyshev series coefficients, ordered
  314. from lowest order term to highest.
  315. Returns
  316. -------
  317. pol : ndarray
  318. 1-D array containing the coefficients of the equivalent polynomial
  319. (relative to the "standard" basis) ordered from lowest order term
  320. to highest.
  321. See Also
  322. --------
  323. poly2cheb
  324. Notes
  325. -----
  326. The easy way to do conversions between polynomial basis sets
  327. is to use the convert method of a class instance.
  328. Examples
  329. --------
  330. >>> from numpy import polynomial as P
  331. >>> c = P.Chebyshev(range(4))
  332. >>> c
  333. Chebyshev([ 0., 1., 2., 3.], [-1., 1.])
  334. >>> p = c.convert(kind=P.Polynomial)
  335. >>> p
  336. Polynomial([ -2., -8., 4., 12.], [-1., 1.])
  337. >>> P.cheb2poly(range(4))
  338. array([ -2., -8., 4., 12.])
  339. """
  340. from .polynomial import polyadd, polysub, polymulx
  341. [c] = pu.as_series([c])
  342. n = len(c)
  343. if n < 3:
  344. return c
  345. else:
  346. c0 = c[-2]
  347. c1 = c[-1]
  348. # i is the current degree of c1
  349. for i in range(n - 1, 1, -1):
  350. tmp = c0
  351. c0 = polysub(c[i - 2], c1)
  352. c1 = polyadd(tmp, polymulx(c1)*2)
  353. return polyadd(c0, polymulx(c1))
  354. #
  355. # These are constant arrays are of integer type so as to be compatible
  356. # with the widest range of other types, such as Decimal.
  357. #
  358. # Chebyshev default domain.
  359. chebdomain = np.array([-1, 1])
  360. # Chebyshev coefficients representing zero.
  361. chebzero = np.array([0])
  362. # Chebyshev coefficients representing one.
  363. chebone = np.array([1])
  364. # Chebyshev coefficients representing the identity x.
  365. chebx = np.array([0, 1])
  366. def chebline(off, scl):
  367. """
  368. Chebyshev series whose graph is a straight line.
  369. Parameters
  370. ----------
  371. off, scl : scalars
  372. The specified line is given by ``off + scl*x``.
  373. Returns
  374. -------
  375. y : ndarray
  376. This module's representation of the Chebyshev series for
  377. ``off + scl*x``.
  378. See Also
  379. --------
  380. polyline
  381. Examples
  382. --------
  383. >>> import numpy.polynomial.chebyshev as C
  384. >>> C.chebline(3,2)
  385. array([3, 2])
  386. >>> C.chebval(-3, C.chebline(3,2)) # should be -3
  387. -3.0
  388. """
  389. if scl != 0:
  390. return np.array([off, scl])
  391. else:
  392. return np.array([off])
  393. def chebfromroots(roots):
  394. """
  395. Generate a Chebyshev series with given roots.
  396. The function returns the coefficients of the polynomial
  397. .. math:: p(x) = (x - r_0) * (x - r_1) * ... * (x - r_n),
  398. in Chebyshev form, where the `r_n` are the roots specified in `roots`.
  399. If a zero has multiplicity n, then it must appear in `roots` n times.
  400. For instance, if 2 is a root of multiplicity three and 3 is a root of
  401. multiplicity 2, then `roots` looks something like [2, 2, 2, 3, 3]. The
  402. roots can appear in any order.
  403. If the returned coefficients are `c`, then
  404. .. math:: p(x) = c_0 + c_1 * T_1(x) + ... + c_n * T_n(x)
  405. The coefficient of the last term is not generally 1 for monic
  406. polynomials in Chebyshev form.
  407. Parameters
  408. ----------
  409. roots : array_like
  410. Sequence containing the roots.
  411. Returns
  412. -------
  413. out : ndarray
  414. 1-D array of coefficients. If all roots are real then `out` is a
  415. real array, if some of the roots are complex, then `out` is complex
  416. even if all the coefficients in the result are real (see Examples
  417. below).
  418. See Also
  419. --------
  420. polyfromroots, legfromroots, lagfromroots, hermfromroots,
  421. hermefromroots.
  422. Examples
  423. --------
  424. >>> import numpy.polynomial.chebyshev as C
  425. >>> C.chebfromroots((-1,0,1)) # x^3 - x relative to the standard basis
  426. array([ 0. , -0.25, 0. , 0.25])
  427. >>> j = complex(0,1)
  428. >>> C.chebfromroots((-j,j)) # x^2 + 1 relative to the standard basis
  429. array([ 1.5+0.j, 0.0+0.j, 0.5+0.j])
  430. """
  431. if len(roots) == 0:
  432. return np.ones(1)
  433. else:
  434. [roots] = pu.as_series([roots], trim=False)
  435. roots.sort()
  436. p = [chebline(-r, 1) for r in roots]
  437. n = len(p)
  438. while n > 1:
  439. m, r = divmod(n, 2)
  440. tmp = [chebmul(p[i], p[i+m]) for i in range(m)]
  441. if r:
  442. tmp[0] = chebmul(tmp[0], p[-1])
  443. p = tmp
  444. n = m
  445. return p[0]
  446. def chebadd(c1, c2):
  447. """
  448. Add one Chebyshev series to another.
  449. Returns the sum of two Chebyshev series `c1` + `c2`. The arguments
  450. are sequences of coefficients ordered from lowest order term to
  451. highest, i.e., [1,2,3] represents the series ``T_0 + 2*T_1 + 3*T_2``.
  452. Parameters
  453. ----------
  454. c1, c2 : array_like
  455. 1-D arrays of Chebyshev series coefficients ordered from low to
  456. high.
  457. Returns
  458. -------
  459. out : ndarray
  460. Array representing the Chebyshev series of their sum.
  461. See Also
  462. --------
  463. chebsub, chebmul, chebdiv, chebpow
  464. Notes
  465. -----
  466. Unlike multiplication, division, etc., the sum of two Chebyshev series
  467. is a Chebyshev series (without having to "reproject" the result onto
  468. the basis set) so addition, just like that of "standard" polynomials,
  469. is simply "component-wise."
  470. Examples
  471. --------
  472. >>> from numpy.polynomial import chebyshev as C
  473. >>> c1 = (1,2,3)
  474. >>> c2 = (3,2,1)
  475. >>> C.chebadd(c1,c2)
  476. array([ 4., 4., 4.])
  477. """
  478. # c1, c2 are trimmed copies
  479. [c1, c2] = pu.as_series([c1, c2])
  480. if len(c1) > len(c2):
  481. c1[:c2.size] += c2
  482. ret = c1
  483. else:
  484. c2[:c1.size] += c1
  485. ret = c2
  486. return pu.trimseq(ret)
  487. def chebsub(c1, c2):
  488. """
  489. Subtract one Chebyshev series from another.
  490. Returns the difference of two Chebyshev series `c1` - `c2`. The
  491. sequences of coefficients are from lowest order term to highest, i.e.,
  492. [1,2,3] represents the series ``T_0 + 2*T_1 + 3*T_2``.
  493. Parameters
  494. ----------
  495. c1, c2 : array_like
  496. 1-D arrays of Chebyshev series coefficients ordered from low to
  497. high.
  498. Returns
  499. -------
  500. out : ndarray
  501. Of Chebyshev series coefficients representing their difference.
  502. See Also
  503. --------
  504. chebadd, chebmul, chebdiv, chebpow
  505. Notes
  506. -----
  507. Unlike multiplication, division, etc., the difference of two Chebyshev
  508. series is a Chebyshev series (without having to "reproject" the result
  509. onto the basis set) so subtraction, just like that of "standard"
  510. polynomials, is simply "component-wise."
  511. Examples
  512. --------
  513. >>> from numpy.polynomial import chebyshev as C
  514. >>> c1 = (1,2,3)
  515. >>> c2 = (3,2,1)
  516. >>> C.chebsub(c1,c2)
  517. array([-2., 0., 2.])
  518. >>> C.chebsub(c2,c1) # -C.chebsub(c1,c2)
  519. array([ 2., 0., -2.])
  520. """
  521. # c1, c2 are trimmed copies
  522. [c1, c2] = pu.as_series([c1, c2])
  523. if len(c1) > len(c2):
  524. c1[:c2.size] -= c2
  525. ret = c1
  526. else:
  527. c2 = -c2
  528. c2[:c1.size] += c1
  529. ret = c2
  530. return pu.trimseq(ret)
  531. def chebmulx(c):
  532. """Multiply a Chebyshev series by x.
  533. Multiply the polynomial `c` by x, where x is the independent
  534. variable.
  535. Parameters
  536. ----------
  537. c : array_like
  538. 1-D array of Chebyshev series coefficients ordered from low to
  539. high.
  540. Returns
  541. -------
  542. out : ndarray
  543. Array representing the result of the multiplication.
  544. Notes
  545. -----
  546. .. versionadded:: 1.5.0
  547. """
  548. # c is a trimmed copy
  549. [c] = pu.as_series([c])
  550. # The zero series needs special treatment
  551. if len(c) == 1 and c[0] == 0:
  552. return c
  553. prd = np.empty(len(c) + 1, dtype=c.dtype)
  554. prd[0] = c[0]*0
  555. prd[1] = c[0]
  556. if len(c) > 1:
  557. tmp = c[1:]/2
  558. prd[2:] = tmp
  559. prd[0:-2] += tmp
  560. return prd
  561. def chebmul(c1, c2):
  562. """
  563. Multiply one Chebyshev series by another.
  564. Returns the product of two Chebyshev series `c1` * `c2`. The arguments
  565. are sequences of coefficients, from lowest order "term" to highest,
  566. e.g., [1,2,3] represents the series ``T_0 + 2*T_1 + 3*T_2``.
  567. Parameters
  568. ----------
  569. c1, c2 : array_like
  570. 1-D arrays of Chebyshev series coefficients ordered from low to
  571. high.
  572. Returns
  573. -------
  574. out : ndarray
  575. Of Chebyshev series coefficients representing their product.
  576. See Also
  577. --------
  578. chebadd, chebsub, chebdiv, chebpow
  579. Notes
  580. -----
  581. In general, the (polynomial) product of two C-series results in terms
  582. that are not in the Chebyshev polynomial basis set. Thus, to express
  583. the product as a C-series, it is typically necessary to "reproject"
  584. the product onto said basis set, which typically produces
  585. "unintuitive live" (but correct) results; see Examples section below.
  586. Examples
  587. --------
  588. >>> from numpy.polynomial import chebyshev as C
  589. >>> c1 = (1,2,3)
  590. >>> c2 = (3,2,1)
  591. >>> C.chebmul(c1,c2) # multiplication requires "reprojection"
  592. array([ 6.5, 12. , 12. , 4. , 1.5])
  593. """
  594. # c1, c2 are trimmed copies
  595. [c1, c2] = pu.as_series([c1, c2])
  596. z1 = _cseries_to_zseries(c1)
  597. z2 = _cseries_to_zseries(c2)
  598. prd = _zseries_mul(z1, z2)
  599. ret = _zseries_to_cseries(prd)
  600. return pu.trimseq(ret)
  601. def chebdiv(c1, c2):
  602. """
  603. Divide one Chebyshev series by another.
  604. Returns the quotient-with-remainder of two Chebyshev series
  605. `c1` / `c2`. The arguments are sequences of coefficients from lowest
  606. order "term" to highest, e.g., [1,2,3] represents the series
  607. ``T_0 + 2*T_1 + 3*T_2``.
  608. Parameters
  609. ----------
  610. c1, c2 : array_like
  611. 1-D arrays of Chebyshev series coefficients ordered from low to
  612. high.
  613. Returns
  614. -------
  615. [quo, rem] : ndarrays
  616. Of Chebyshev series coefficients representing the quotient and
  617. remainder.
  618. See Also
  619. --------
  620. chebadd, chebsub, chebmul, chebpow
  621. Notes
  622. -----
  623. In general, the (polynomial) division of one C-series by another
  624. results in quotient and remainder terms that are not in the Chebyshev
  625. polynomial basis set. Thus, to express these results as C-series, it
  626. is typically necessary to "reproject" the results onto said basis
  627. set, which typically produces "unintuitive" (but correct) results;
  628. see Examples section below.
  629. Examples
  630. --------
  631. >>> from numpy.polynomial import chebyshev as C
  632. >>> c1 = (1,2,3)
  633. >>> c2 = (3,2,1)
  634. >>> C.chebdiv(c1,c2) # quotient "intuitive," remainder not
  635. (array([ 3.]), array([-8., -4.]))
  636. >>> c2 = (0,1,2,3)
  637. >>> C.chebdiv(c2,c1) # neither "intuitive"
  638. (array([ 0., 2.]), array([-2., -4.]))
  639. """
  640. # c1, c2 are trimmed copies
  641. [c1, c2] = pu.as_series([c1, c2])
  642. if c2[-1] == 0:
  643. raise ZeroDivisionError()
  644. lc1 = len(c1)
  645. lc2 = len(c2)
  646. if lc1 < lc2:
  647. return c1[:1]*0, c1
  648. elif lc2 == 1:
  649. return c1/c2[-1], c1[:1]*0
  650. else:
  651. z1 = _cseries_to_zseries(c1)
  652. z2 = _cseries_to_zseries(c2)
  653. quo, rem = _zseries_div(z1, z2)
  654. quo = pu.trimseq(_zseries_to_cseries(quo))
  655. rem = pu.trimseq(_zseries_to_cseries(rem))
  656. return quo, rem
  657. def chebpow(c, pow, maxpower=16):
  658. """Raise a Chebyshev series to a power.
  659. Returns the Chebyshev series `c` raised to the power `pow`. The
  660. argument `c` is a sequence of coefficients ordered from low to high.
  661. i.e., [1,2,3] is the series ``T_0 + 2*T_1 + 3*T_2.``
  662. Parameters
  663. ----------
  664. c : array_like
  665. 1-D array of Chebyshev series coefficients ordered from low to
  666. high.
  667. pow : integer
  668. Power to which the series will be raised
  669. maxpower : integer, optional
  670. Maximum power allowed. This is mainly to limit growth of the series
  671. to unmanageable size. Default is 16
  672. Returns
  673. -------
  674. coef : ndarray
  675. Chebyshev series of power.
  676. See Also
  677. --------
  678. chebadd, chebsub, chebmul, chebdiv
  679. Examples
  680. --------
  681. """
  682. # c is a trimmed copy
  683. [c] = pu.as_series([c])
  684. power = int(pow)
  685. if power != pow or power < 0:
  686. raise ValueError("Power must be a non-negative integer.")
  687. elif maxpower is not None and power > maxpower:
  688. raise ValueError("Power is too large")
  689. elif power == 0:
  690. return np.array([1], dtype=c.dtype)
  691. elif power == 1:
  692. return c
  693. else:
  694. # This can be made more efficient by using powers of two
  695. # in the usual way.
  696. zs = _cseries_to_zseries(c)
  697. prd = zs
  698. for i in range(2, power + 1):
  699. prd = np.convolve(prd, zs)
  700. return _zseries_to_cseries(prd)
  701. def chebder(c, m=1, scl=1, axis=0):
  702. """
  703. Differentiate a Chebyshev series.
  704. Returns the Chebyshev series coefficients `c` differentiated `m` times
  705. along `axis`. At each iteration the result is multiplied by `scl` (the
  706. scaling factor is for use in a linear change of variable). The argument
  707. `c` is an array of coefficients from low to high degree along each
  708. axis, e.g., [1,2,3] represents the series ``1*T_0 + 2*T_1 + 3*T_2``
  709. while [[1,2],[1,2]] represents ``1*T_0(x)*T_0(y) + 1*T_1(x)*T_0(y) +
  710. 2*T_0(x)*T_1(y) + 2*T_1(x)*T_1(y)`` if axis=0 is ``x`` and axis=1 is
  711. ``y``.
  712. Parameters
  713. ----------
  714. c : array_like
  715. Array of Chebyshev series coefficients. If c is multidimensional
  716. the different axis correspond to different variables with the
  717. degree in each axis given by the corresponding index.
  718. m : int, optional
  719. Number of derivatives taken, must be non-negative. (Default: 1)
  720. scl : scalar, optional
  721. Each differentiation is multiplied by `scl`. The end result is
  722. multiplication by ``scl**m``. This is for use in a linear change of
  723. variable. (Default: 1)
  724. axis : int, optional
  725. Axis over which the derivative is taken. (Default: 0).
  726. .. versionadded:: 1.7.0
  727. Returns
  728. -------
  729. der : ndarray
  730. Chebyshev series of the derivative.
  731. See Also
  732. --------
  733. chebint
  734. Notes
  735. -----
  736. In general, the result of differentiating a C-series needs to be
  737. "reprojected" onto the C-series basis set. Thus, typically, the
  738. result of this function is "unintuitive," albeit correct; see Examples
  739. section below.
  740. Examples
  741. --------
  742. >>> from numpy.polynomial import chebyshev as C
  743. >>> c = (1,2,3,4)
  744. >>> C.chebder(c)
  745. array([ 14., 12., 24.])
  746. >>> C.chebder(c,3)
  747. array([ 96.])
  748. >>> C.chebder(c,scl=-1)
  749. array([-14., -12., -24.])
  750. >>> C.chebder(c,2,-1)
  751. array([ 12., 96.])
  752. """
  753. c = np.array(c, ndmin=1, copy=1)
  754. if c.dtype.char in '?bBhHiIlLqQpP':
  755. c = c.astype(np.double)
  756. cnt, iaxis = [int(t) for t in [m, axis]]
  757. if cnt != m:
  758. raise ValueError("The order of derivation must be integer")
  759. if cnt < 0:
  760. raise ValueError("The order of derivation must be non-negative")
  761. if iaxis != axis:
  762. raise ValueError("The axis must be integer")
  763. if not -c.ndim <= iaxis < c.ndim:
  764. raise ValueError("The axis is out of range")
  765. if iaxis < 0:
  766. iaxis += c.ndim
  767. if cnt == 0:
  768. return c
  769. c = np.rollaxis(c, iaxis)
  770. n = len(c)
  771. if cnt >= n:
  772. c = c[:1]*0
  773. else:
  774. for i in range(cnt):
  775. n = n - 1
  776. c *= scl
  777. der = np.empty((n,) + c.shape[1:], dtype=c.dtype)
  778. for j in range(n, 2, -1):
  779. der[j - 1] = (2*j)*c[j]
  780. c[j - 2] += (j*c[j])/(j - 2)
  781. if n > 1:
  782. der[1] = 4*c[2]
  783. der[0] = c[1]
  784. c = der
  785. c = np.rollaxis(c, 0, iaxis + 1)
  786. return c
  787. def chebint(c, m=1, k=[], lbnd=0, scl=1, axis=0):
  788. """
  789. Integrate a Chebyshev series.
  790. Returns the Chebyshev series coefficients `c` integrated `m` times from
  791. `lbnd` along `axis`. At each iteration the resulting series is
  792. **multiplied** by `scl` and an integration constant, `k`, is added.
  793. The scaling factor is for use in a linear change of variable. ("Buyer
  794. beware": note that, depending on what one is doing, one may want `scl`
  795. to be the reciprocal of what one might expect; for more information,
  796. see the Notes section below.) The argument `c` is an array of
  797. coefficients from low to high degree along each axis, e.g., [1,2,3]
  798. represents the series ``T_0 + 2*T_1 + 3*T_2`` while [[1,2],[1,2]]
  799. represents ``1*T_0(x)*T_0(y) + 1*T_1(x)*T_0(y) + 2*T_0(x)*T_1(y) +
  800. 2*T_1(x)*T_1(y)`` if axis=0 is ``x`` and axis=1 is ``y``.
  801. Parameters
  802. ----------
  803. c : array_like
  804. Array of Chebyshev series coefficients. If c is multidimensional
  805. the different axis correspond to different variables with the
  806. degree in each axis given by the corresponding index.
  807. m : int, optional
  808. Order of integration, must be positive. (Default: 1)
  809. k : {[], list, scalar}, optional
  810. Integration constant(s). The value of the first integral at zero
  811. is the first value in the list, the value of the second integral
  812. at zero is the second value, etc. If ``k == []`` (the default),
  813. all constants are set to zero. If ``m == 1``, a single scalar can
  814. be given instead of a list.
  815. lbnd : scalar, optional
  816. The lower bound of the integral. (Default: 0)
  817. scl : scalar, optional
  818. Following each integration the result is *multiplied* by `scl`
  819. before the integration constant is added. (Default: 1)
  820. axis : int, optional
  821. Axis over which the integral is taken. (Default: 0).
  822. .. versionadded:: 1.7.0
  823. Returns
  824. -------
  825. S : ndarray
  826. C-series coefficients of the integral.
  827. Raises
  828. ------
  829. ValueError
  830. If ``m < 1``, ``len(k) > m``, ``np.isscalar(lbnd) == False``, or
  831. ``np.isscalar(scl) == False``.
  832. See Also
  833. --------
  834. chebder
  835. Notes
  836. -----
  837. Note that the result of each integration is *multiplied* by `scl`.
  838. Why is this important to note? Say one is making a linear change of
  839. variable :math:`u = ax + b` in an integral relative to `x`. Then
  840. .. math::`dx = du/a`, so one will need to set `scl` equal to
  841. :math:`1/a`- perhaps not what one would have first thought.
  842. Also note that, in general, the result of integrating a C-series needs
  843. to be "reprojected" onto the C-series basis set. Thus, typically,
  844. the result of this function is "unintuitive," albeit correct; see
  845. Examples section below.
  846. Examples
  847. --------
  848. >>> from numpy.polynomial import chebyshev as C
  849. >>> c = (1,2,3)
  850. >>> C.chebint(c)
  851. array([ 0.5, -0.5, 0.5, 0.5])
  852. >>> C.chebint(c,3)
  853. array([ 0.03125 , -0.1875 , 0.04166667, -0.05208333, 0.01041667,
  854. 0.00625 ])
  855. >>> C.chebint(c, k=3)
  856. array([ 3.5, -0.5, 0.5, 0.5])
  857. >>> C.chebint(c,lbnd=-2)
  858. array([ 8.5, -0.5, 0.5, 0.5])
  859. >>> C.chebint(c,scl=-2)
  860. array([-1., 1., -1., -1.])
  861. """
  862. c = np.array(c, ndmin=1, copy=1)
  863. if c.dtype.char in '?bBhHiIlLqQpP':
  864. c = c.astype(np.double)
  865. if not np.iterable(k):
  866. k = [k]
  867. cnt, iaxis = [int(t) for t in [m, axis]]
  868. if cnt != m:
  869. raise ValueError("The order of integration must be integer")
  870. if cnt < 0:
  871. raise ValueError("The order of integration must be non-negative")
  872. if len(k) > cnt:
  873. raise ValueError("Too many integration constants")
  874. if iaxis != axis:
  875. raise ValueError("The axis must be integer")
  876. if not -c.ndim <= iaxis < c.ndim:
  877. raise ValueError("The axis is out of range")
  878. if iaxis < 0:
  879. iaxis += c.ndim
  880. if cnt == 0:
  881. return c
  882. c = np.rollaxis(c, iaxis)
  883. k = list(k) + [0]*(cnt - len(k))
  884. for i in range(cnt):
  885. n = len(c)
  886. c *= scl
  887. if n == 1 and np.all(c[0] == 0):
  888. c[0] += k[i]
  889. else:
  890. tmp = np.empty((n + 1,) + c.shape[1:], dtype=c.dtype)
  891. tmp[0] = c[0]*0
  892. tmp[1] = c[0]
  893. if n > 1:
  894. tmp[2] = c[1]/4
  895. for j in range(2, n):
  896. t = c[j]/(2*j + 1)
  897. tmp[j + 1] = c[j]/(2*(j + 1))
  898. tmp[j - 1] -= c[j]/(2*(j - 1))
  899. tmp[0] += k[i] - chebval(lbnd, tmp)
  900. c = tmp
  901. c = np.rollaxis(c, 0, iaxis + 1)
  902. return c
  903. def chebval(x, c, tensor=True):
  904. """
  905. Evaluate a Chebyshev series at points x.
  906. If `c` is of length `n + 1`, this function returns the value:
  907. .. math:: p(x) = c_0 * T_0(x) + c_1 * T_1(x) + ... + c_n * T_n(x)
  908. The parameter `x` is converted to an array only if it is a tuple or a
  909. list, otherwise it is treated as a scalar. In either case, either `x`
  910. or its elements must support multiplication and addition both with
  911. themselves and with the elements of `c`.
  912. If `c` is a 1-D array, then `p(x)` will have the same shape as `x`. If
  913. `c` is multidimensional, then the shape of the result depends on the
  914. value of `tensor`. If `tensor` is true the shape will be c.shape[1:] +
  915. x.shape. If `tensor` is false the shape will be c.shape[1:]. Note that
  916. scalars have shape (,).
  917. Trailing zeros in the coefficients will be used in the evaluation, so
  918. they should be avoided if efficiency is a concern.
  919. Parameters
  920. ----------
  921. x : array_like, compatible object
  922. If `x` is a list or tuple, it is converted to an ndarray, otherwise
  923. it is left unchanged and treated as a scalar. In either case, `x`
  924. or its elements must support addition and multiplication with
  925. with themselves and with the elements of `c`.
  926. c : array_like
  927. Array of coefficients ordered so that the coefficients for terms of
  928. degree n are contained in c[n]. If `c` is multidimensional the
  929. remaining indices enumerate multiple polynomials. In the two
  930. dimensional case the coefficients may be thought of as stored in
  931. the columns of `c`.
  932. tensor : boolean, optional
  933. If True, the shape of the coefficient array is extended with ones
  934. on the right, one for each dimension of `x`. Scalars have dimension 0
  935. for this action. The result is that every column of coefficients in
  936. `c` is evaluated for every element of `x`. If False, `x` is broadcast
  937. over the columns of `c` for the evaluation. This keyword is useful
  938. when `c` is multidimensional. The default value is True.
  939. .. versionadded:: 1.7.0
  940. Returns
  941. -------
  942. values : ndarray, algebra_like
  943. The shape of the return value is described above.
  944. See Also
  945. --------
  946. chebval2d, chebgrid2d, chebval3d, chebgrid3d
  947. Notes
  948. -----
  949. The evaluation uses Clenshaw recursion, aka synthetic division.
  950. Examples
  951. --------
  952. """
  953. c = np.array(c, ndmin=1, copy=1)
  954. if c.dtype.char in '?bBhHiIlLqQpP':
  955. c = c.astype(np.double)
  956. if isinstance(x, (tuple, list)):
  957. x = np.asarray(x)
  958. if isinstance(x, np.ndarray) and tensor:
  959. c = c.reshape(c.shape + (1,)*x.ndim)
  960. if len(c) == 1:
  961. c0 = c[0]
  962. c1 = 0
  963. elif len(c) == 2:
  964. c0 = c[0]
  965. c1 = c[1]
  966. else:
  967. x2 = 2*x
  968. c0 = c[-2]
  969. c1 = c[-1]
  970. for i in range(3, len(c) + 1):
  971. tmp = c0
  972. c0 = c[-i] - c1
  973. c1 = tmp + c1*x2
  974. return c0 + c1*x
  975. def chebval2d(x, y, c):
  976. """
  977. Evaluate a 2-D Chebyshev series at points (x, y).
  978. This function returns the values:
  979. .. math:: p(x,y) = \\sum_{i,j} c_{i,j} * T_i(x) * T_j(y)
  980. The parameters `x` and `y` are converted to arrays only if they are
  981. tuples or a lists, otherwise they are treated as a scalars and they
  982. must have the same shape after conversion. In either case, either `x`
  983. and `y` or their elements must support multiplication and addition both
  984. with themselves and with the elements of `c`.
  985. If `c` is a 1-D array a one is implicitly appended to its shape to make
  986. it 2-D. The shape of the result will be c.shape[2:] + x.shape.
  987. Parameters
  988. ----------
  989. x, y : array_like, compatible objects
  990. The two dimensional series is evaluated at the points `(x, y)`,
  991. where `x` and `y` must have the same shape. If `x` or `y` is a list
  992. or tuple, it is first converted to an ndarray, otherwise it is left
  993. unchanged and if it isn't an ndarray it is treated as a scalar.
  994. c : array_like
  995. Array of coefficients ordered so that the coefficient of the term
  996. of multi-degree i,j is contained in ``c[i,j]``. If `c` has
  997. dimension greater than 2 the remaining indices enumerate multiple
  998. sets of coefficients.
  999. Returns
  1000. -------
  1001. values : ndarray, compatible object
  1002. The values of the two dimensional Chebyshev series at points formed
  1003. from pairs of corresponding values from `x` and `y`.
  1004. See Also
  1005. --------
  1006. chebval, chebgrid2d, chebval3d, chebgrid3d
  1007. Notes
  1008. -----
  1009. .. versionadded::1.7.0
  1010. """
  1011. try:
  1012. x, y = np.array((x, y), copy=0)
  1013. except:
  1014. raise ValueError('x, y are incompatible')
  1015. c = chebval(x, c)
  1016. c = chebval(y, c, tensor=False)
  1017. return c
  1018. def chebgrid2d(x, y, c):
  1019. """
  1020. Evaluate a 2-D Chebyshev series on the Cartesian product of x and y.
  1021. This function returns the values:
  1022. .. math:: p(a,b) = \sum_{i,j} c_{i,j} * T_i(a) * T_j(b),
  1023. where the points `(a, b)` consist of all pairs formed by taking
  1024. `a` from `x` and `b` from `y`. The resulting points form a grid with
  1025. `x` in the first dimension and `y` in the second.
  1026. The parameters `x` and `y` are converted to arrays only if they are
  1027. tuples or a lists, otherwise they are treated as a scalars. In either
  1028. case, either `x` and `y` or their elements must support multiplication
  1029. and addition both with themselves and with the elements of `c`.
  1030. If `c` has fewer than two dimensions, ones are implicitly appended to
  1031. its shape to make it 2-D. The shape of the result will be c.shape[2:] +
  1032. x.shape + y.shape.
  1033. Parameters
  1034. ----------
  1035. x, y : array_like, compatible objects
  1036. The two dimensional series is evaluated at the points in the
  1037. Cartesian product of `x` and `y`. If `x` or `y` is a list or
  1038. tuple, it is first converted to an ndarray, otherwise it is left
  1039. unchanged and, if it isn't an ndarray, it is treated as a scalar.
  1040. c : array_like
  1041. Array of coefficients ordered so that the coefficient of the term of
  1042. multi-degree i,j is contained in `c[i,j]`. If `c` has dimension
  1043. greater than two the remaining indices enumerate multiple sets of
  1044. coefficients.
  1045. Returns
  1046. -------
  1047. values : ndarray, compatible object
  1048. The values of the two dimensional Chebyshev series at points in the
  1049. Cartesian product of `x` and `y`.
  1050. See Also
  1051. --------
  1052. chebval, chebval2d, chebval3d, chebgrid3d
  1053. Notes
  1054. -----
  1055. .. versionadded::1.7.0
  1056. """
  1057. c = chebval(x, c)
  1058. c = chebval(y, c)
  1059. return c
  1060. def chebval3d(x, y, z, c):
  1061. """
  1062. Evaluate a 3-D Chebyshev series at points (x, y, z).
  1063. This function returns the values:
  1064. .. math:: p(x,y,z) = \\sum_{i,j,k} c_{i,j,k} * T_i(x) * T_j(y) * T_k(z)
  1065. The parameters `x`, `y`, and `z` are converted to arrays only if
  1066. they are tuples or a lists, otherwise they are treated as a scalars and
  1067. they must have the same shape after conversion. In either case, either
  1068. `x`, `y`, and `z` or their elements must support multiplication and
  1069. addition both with themselves and with the elements of `c`.
  1070. If `c` has fewer than 3 dimensions, ones are implicitly appended to its
  1071. shape to make it 3-D. The shape of the result will be c.shape[3:] +
  1072. x.shape.
  1073. Parameters
  1074. ----------
  1075. x, y, z : array_like, compatible object
  1076. The three dimensional series is evaluated at the points
  1077. `(x, y, z)`, where `x`, `y`, and `z` must have the same shape. If
  1078. any of `x`, `y`, or `z` is a list or tuple, it is first converted
  1079. to an ndarray, otherwise it is left unchanged and if it isn't an
  1080. ndarray it is treated as a scalar.
  1081. c : array_like
  1082. Array of coefficients ordered so that the coefficient of the term of
  1083. multi-degree i,j,k is contained in ``c[i,j,k]``. If `c` has dimension
  1084. greater than 3 the remaining indices enumerate multiple sets of
  1085. coefficients.
  1086. Returns
  1087. -------
  1088. values : ndarray, compatible object
  1089. The values of the multidimensional polynomial on points formed with
  1090. triples of corresponding values from `x`, `y`, and `z`.
  1091. See Also
  1092. --------
  1093. chebval, chebval2d, chebgrid2d, chebgrid3d
  1094. Notes
  1095. -----
  1096. .. versionadded::1.7.0
  1097. """
  1098. try:
  1099. x, y, z = np.array((x, y, z), copy=0)
  1100. except:
  1101. raise ValueError('x, y, z are incompatible')
  1102. c = chebval(x, c)
  1103. c = chebval(y, c, tensor=False)
  1104. c = chebval(z, c, tensor=False)
  1105. return c
  1106. def chebgrid3d(x, y, z, c):
  1107. """
  1108. Evaluate a 3-D Chebyshev series on the Cartesian product of x, y, and z.
  1109. This function returns the values:
  1110. .. math:: p(a,b,c) = \\sum_{i,j,k} c_{i,j,k} * T_i(a) * T_j(b) * T_k(c)
  1111. where the points `(a, b, c)` consist of all triples formed by taking
  1112. `a` from `x`, `b` from `y`, and `c` from `z`. The resulting points form
  1113. a grid with `x` in the first dimension, `y` in the second, and `z` in
  1114. the third.
  1115. The parameters `x`, `y`, and `z` are converted to arrays only if they
  1116. are tuples or a lists, otherwise they are treated as a scalars. In
  1117. either case, either `x`, `y`, and `z` or their elements must support
  1118. multiplication and addition both with themselves and with the elements
  1119. of `c`.
  1120. If `c` has fewer than three dimensions, ones are implicitly appended to
  1121. its shape to make it 3-D. The shape of the result will be c.shape[3:] +
  1122. x.shape + y.shape + z.shape.
  1123. Parameters
  1124. ----------
  1125. x, y, z : array_like, compatible objects
  1126. The three dimensional series is evaluated at the points in the
  1127. Cartesian product of `x`, `y`, and `z`. If `x`,`y`, or `z` is a
  1128. list or tuple, it is first converted to an ndarray, otherwise it is
  1129. left unchanged and, if it isn't an ndarray, it is treated as a
  1130. scalar.
  1131. c : array_like
  1132. Array of coefficients ordered so that the coefficients for terms of
  1133. degree i,j are contained in ``c[i,j]``. If `c` has dimension
  1134. greater than two the remaining indices enumerate multiple sets of
  1135. coefficients.
  1136. Returns
  1137. -------
  1138. values : ndarray, compatible object
  1139. The values of the two dimensional polynomial at points in the Cartesian
  1140. product of `x` and `y`.
  1141. See Also
  1142. --------
  1143. chebval, chebval2d, chebgrid2d, chebval3d
  1144. Notes
  1145. -----
  1146. .. versionadded::1.7.0
  1147. """
  1148. c = chebval(x, c)
  1149. c = chebval(y, c)
  1150. c = chebval(z, c)
  1151. return c
  1152. def chebvander(x, deg):
  1153. """Pseudo-Vandermonde matrix of given degree.
  1154. Returns the pseudo-Vandermonde matrix of degree `deg` and sample points
  1155. `x`. The pseudo-Vandermonde matrix is defined by
  1156. .. math:: V[..., i] = T_i(x),
  1157. where `0 <= i <= deg`. The leading indices of `V` index the elements of
  1158. `x` and the last index is the degree of the Chebyshev polynomial.
  1159. If `c` is a 1-D array of coefficients of length `n + 1` and `V` is the
  1160. matrix ``V = chebvander(x, n)``, then ``np.dot(V, c)`` and
  1161. ``chebval(x, c)`` are the same up to roundoff. This equivalence is
  1162. useful both for least squares fitting and for the evaluation of a large
  1163. number of Chebyshev series of the same degree and sample points.
  1164. Parameters
  1165. ----------
  1166. x : array_like
  1167. Array of points. The dtype is converted to float64 or complex128
  1168. depending on whether any of the elements are complex. If `x` is
  1169. scalar it is converted to a 1-D array.
  1170. deg : int
  1171. Degree of the resulting matrix.
  1172. Returns
  1173. -------
  1174. vander : ndarray
  1175. The pseudo Vandermonde matrix. The shape of the returned matrix is
  1176. ``x.shape + (deg + 1,)``, where The last index is the degree of the
  1177. corresponding Chebyshev polynomial. The dtype will be the same as
  1178. the converted `x`.
  1179. """
  1180. ideg = int(deg)
  1181. if ideg != deg:
  1182. raise ValueError("deg must be integer")
  1183. if ideg < 0:
  1184. raise ValueError("deg must be non-negative")
  1185. x = np.array(x, copy=0, ndmin=1) + 0.0
  1186. dims = (ideg + 1,) + x.shape
  1187. dtyp = x.dtype
  1188. v = np.empty(dims, dtype=dtyp)
  1189. # Use forward recursion to generate the entries.
  1190. v[0] = x*0 + 1
  1191. if ideg > 0:
  1192. x2 = 2*x
  1193. v[1] = x
  1194. for i in range(2, ideg + 1):
  1195. v[i] = v[i-1]*x2 - v[i-2]
  1196. return np.rollaxis(v, 0, v.ndim)
  1197. def chebvander2d(x, y, deg):
  1198. """Pseudo-Vandermonde matrix of given degrees.
  1199. Returns the pseudo-Vandermonde matrix of degrees `deg` and sample
  1200. points `(x, y)`. The pseudo-Vandermonde matrix is defined by
  1201. .. math:: V[..., deg[1]*i + j] = T_i(x) * T_j(y),
  1202. where `0 <= i <= deg[0]` and `0 <= j <= deg[1]`. The leading indices of
  1203. `V` index the points `(x, y)` and the last index encodes the degrees of
  1204. the Chebyshev polynomials.
  1205. If ``V = chebvander2d(x, y, [xdeg, ydeg])``, then the columns of `V`
  1206. correspond to the elements of a 2-D coefficient array `c` of shape
  1207. (xdeg + 1, ydeg + 1) in the order
  1208. .. math:: c_{00}, c_{01}, c_{02} ... , c_{10}, c_{11}, c_{12} ...
  1209. and ``np.dot(V, c.flat)`` and ``chebval2d(x, y, c)`` will be the same
  1210. up to roundoff. This equivalence is useful both for least squares
  1211. fitting and for the evaluation of a large number of 2-D Chebyshev
  1212. series of the same degrees and sample points.
  1213. Parameters
  1214. ----------
  1215. x, y : array_like
  1216. Arrays of point coordinates, all of the same shape. The dtypes
  1217. will be converted to either float64 or complex128 depending on
  1218. whether any of the elements are complex. Scalars are converted to
  1219. 1-D arrays.
  1220. deg : list of ints
  1221. List of maximum degrees of the form [x_deg, y_deg].
  1222. Returns
  1223. -------
  1224. vander2d : ndarray
  1225. The shape of the returned matrix is ``x.shape + (order,)``, where
  1226. :math:`order = (deg[0]+1)*(deg([1]+1)`. The dtype will be the same
  1227. as the converted `x` and `y`.
  1228. See Also
  1229. --------
  1230. chebvander, chebvander3d. chebval2d, chebval3d
  1231. Notes
  1232. -----
  1233. .. versionadded::1.7.0
  1234. """
  1235. ideg = [int(d) for d in deg]
  1236. is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
  1237. if is_valid != [1, 1]:
  1238. raise ValueError("degrees must be non-negative integers")
  1239. degx, degy = ideg
  1240. x, y = np.array((x, y), copy=0) + 0.0
  1241. vx = chebvander(x, degx)
  1242. vy = chebvander(y, degy)
  1243. v = vx[..., None]*vy[..., None,:]
  1244. return v.reshape(v.shape[:-2] + (-1,))
  1245. def chebvander3d(x, y, z, deg):
  1246. """Pseudo-Vandermonde matrix of given degrees.
  1247. Returns the pseudo-Vandermonde matrix of degrees `deg` and sample
  1248. points `(x, y, z)`. If `l, m, n` are the given degrees in `x, y, z`,
  1249. then The pseudo-Vandermonde matrix is defined by
  1250. .. math:: V[..., (m+1)(n+1)i + (n+1)j + k] = T_i(x)*T_j(y)*T_k(z),
  1251. where `0 <= i <= l`, `0 <= j <= m`, and `0 <= j <= n`. The leading
  1252. indices of `V` index the points `(x, y, z)` and the last index encodes
  1253. the degrees of the Chebyshev polynomials.
  1254. If ``V = chebvander3d(x, y, z, [xdeg, ydeg, zdeg])``, then the columns
  1255. of `V` correspond to the elements of a 3-D coefficient array `c` of
  1256. shape (xdeg + 1, ydeg + 1, zdeg + 1) in the order
  1257. .. math:: c_{000}, c_{001}, c_{002},... , c_{010}, c_{011}, c_{012},...
  1258. and ``np.dot(V, c.flat)`` and ``chebval3d(x, y, z, c)`` will be the
  1259. same up to roundoff. This equivalence is useful both for least squares
  1260. fitting and for the evaluation of a large number of 3-D Chebyshev
  1261. series of the same degrees and sample points.
  1262. Parameters
  1263. ----------
  1264. x, y, z : array_like
  1265. Arrays of point coordinates, all of the same shape. The dtypes will
  1266. be converted to either float64 or complex128 depending on whether
  1267. any of the elements are complex. Scalars are converted to 1-D
  1268. arrays.
  1269. deg : list of ints
  1270. List of maximum degrees of the form [x_deg, y_deg, z_deg].
  1271. Returns
  1272. -------
  1273. vander3d : ndarray
  1274. The shape of the returned matrix is ``x.shape + (order,)``, where
  1275. :math:`order = (deg[0]+1)*(deg([1]+1)*(deg[2]+1)`. The dtype will
  1276. be the same as the converted `x`, `y`, and `z`.
  1277. See Also
  1278. --------
  1279. chebvander, chebvander3d. chebval2d, chebval3d
  1280. Notes
  1281. -----
  1282. .. versionadded::1.7.0
  1283. """
  1284. ideg = [int(d) for d in deg]
  1285. is_valid = [id == d and id >= 0 for id, d in zip(ideg, deg)]
  1286. if is_valid != [1, 1, 1]:
  1287. raise ValueError("degrees must be non-negative integers")
  1288. degx, degy, degz = ideg
  1289. x, y, z = np.array((x, y, z), copy=0) + 0.0
  1290. vx = chebvander(x, degx)
  1291. vy = chebvander(y, degy)
  1292. vz = chebvander(z, degz)
  1293. v = vx[..., None, None]*vy[..., None,:, None]*vz[..., None, None,:]
  1294. return v.reshape(v.shape[:-3] + (-1,))
  1295. def chebfit(x, y, deg, rcond=None, full=False, w=None):
  1296. """
  1297. Least squares fit of Chebyshev series to data.
  1298. Return the coefficients of a Legendre series of degree `deg` that is the
  1299. least squares fit to the data values `y` given at points `x`. If `y` is
  1300. 1-D the returned coefficients will also be 1-D. If `y` is 2-D multiple
  1301. fits are done, one for each column of `y`, and the resulting
  1302. coefficients are stored in the corresponding columns of a 2-D return.
  1303. The fitted polynomial(s) are in the form
  1304. .. math:: p(x) = c_0 + c_1 * T_1(x) + ... + c_n * T_n(x),
  1305. where `n` is `deg`.
  1306. Parameters
  1307. ----------
  1308. x : array_like, shape (M,)
  1309. x-coordinates of the M sample points ``(x[i], y[i])``.
  1310. y : array_like, shape (M,) or (M, K)
  1311. y-coordinates of the sample points. Several data sets of sample
  1312. points sharing the same x-coordinates can be fitted at once by
  1313. passing in a 2D-array that contains one dataset per column.
  1314. deg : int or 1-D array_like
  1315. Degree(s) of the fitting polynomials. If `deg` is a single integer
  1316. all terms up to and including the `deg`'th term are included in the
  1317. fit. For Numpy versions >= 1.11 a list of integers specifying the
  1318. degrees of the terms to include may be used instead.
  1319. rcond : float, optional
  1320. Relative condition number of the fit. Singular values smaller than
  1321. this relative to the largest singular value will be ignored. The
  1322. default value is len(x)*eps, where eps is the relative precision of
  1323. the float type, about 2e-16 in most cases.
  1324. full : bool, optional
  1325. Switch determining nature of return value. When it is False (the
  1326. default) just the coefficients are returned, when True diagnostic
  1327. information from the singular value decomposition is also returned.
  1328. w : array_like, shape (`M`,), optional
  1329. Weights. If not None, the contribution of each point
  1330. ``(x[i],y[i])`` to the fit is weighted by `w[i]`. Ideally the
  1331. weights are chosen so that the errors of the products ``w[i]*y[i]``
  1332. all have the same variance. The default value is None.
  1333. .. versionadded:: 1.5.0
  1334. Returns
  1335. -------
  1336. coef : ndarray, shape (M,) or (M, K)
  1337. Chebyshev coefficients ordered from low to high. If `y` was 2-D,
  1338. the coefficients for the data in column k of `y` are in column
  1339. `k`.
  1340. [residuals, rank, singular_values, rcond] : list
  1341. These values are only returned if `full` = True
  1342. resid -- sum of squared residuals of the least squares fit
  1343. rank -- the numerical rank of the scaled Vandermonde matrix
  1344. sv -- singular values of the scaled Vandermonde matrix
  1345. rcond -- value of `rcond`.
  1346. For more details, see `linalg.lstsq`.
  1347. Warns
  1348. -----
  1349. RankWarning
  1350. The rank of the coefficient matrix in the least-squares fit is
  1351. deficient. The warning is only raised if `full` = False. The
  1352. warnings can be turned off by
  1353. >>> import warnings
  1354. >>> warnings.simplefilter('ignore', RankWarning)
  1355. See Also
  1356. --------
  1357. polyfit, legfit, lagfit, hermfit, hermefit
  1358. chebval : Evaluates a Chebyshev series.
  1359. chebvander : Vandermonde matrix of Chebyshev series.
  1360. chebweight : Chebyshev weight function.
  1361. linalg.lstsq : Computes a least-squares fit from the matrix.
  1362. scipy.interpolate.UnivariateSpline : Computes spline fits.
  1363. Notes
  1364. -----
  1365. The solution is the coefficients of the Chebyshev series `p` that
  1366. minimizes the sum of the weighted squared errors
  1367. .. math:: E = \\sum_j w_j^2 * |y_j - p(x_j)|^2,
  1368. where :math:`w_j` are the weights. This problem is solved by setting up
  1369. as the (typically) overdetermined matrix equation
  1370. .. math:: V(x) * c = w * y,
  1371. where `V` is the weighted pseudo Vandermonde matrix of `x`, `c` are the
  1372. coefficients to be solved for, `w` are the weights, and `y` are the
  1373. observed values. This equation is then solved using the singular value
  1374. decomposition of `V`.
  1375. If some of the singular values of `V` are so small that they are
  1376. neglected, then a `RankWarning` will be issued. This means that the
  1377. coefficient values may be poorly determined. Using a lower order fit
  1378. will usually get rid of the warning. The `rcond` parameter can also be
  1379. set to a value smaller than its default, but the resulting fit may be
  1380. spurious and have large contributions from roundoff error.
  1381. Fits using Chebyshev series are usually better conditioned than fits
  1382. using power series, but much can depend on the distribution of the
  1383. sample points and the smoothness of the data. If the quality of the fit
  1384. is inadequate splines may be a good alternative.
  1385. References
  1386. ----------
  1387. .. [1] Wikipedia, "Curve fitting",
  1388. http://en.wikipedia.org/wiki/Curve_fitting
  1389. Examples
  1390. --------
  1391. """
  1392. x = np.asarray(x) + 0.0
  1393. y = np.asarray(y) + 0.0
  1394. deg = np.asarray(deg)
  1395. # check arguments.
  1396. if deg.ndim > 1 or deg.dtype.kind not in 'iu' or deg.size == 0:
  1397. raise TypeError("deg must be an int or non-empty 1-D array of int")
  1398. if deg.min() < 0:
  1399. raise ValueError("expected deg >= 0")
  1400. if x.ndim != 1:
  1401. raise TypeError("expected 1D vector for x")
  1402. if x.size == 0:
  1403. raise TypeError("expected non-empty vector for x")
  1404. if y.ndim < 1 or y.ndim > 2:
  1405. raise TypeError("expected 1D or 2D array for y")
  1406. if len(x) != len(y):
  1407. raise TypeError("expected x and y to have same length")
  1408. if deg.ndim == 0:
  1409. lmax = deg
  1410. order = lmax + 1
  1411. van = chebvander(x, lmax)
  1412. else:
  1413. deg = np.sort(deg)
  1414. lmax = deg[-1]
  1415. order = len(deg)
  1416. van = chebvander(x, lmax)[:, deg]
  1417. # set up the least squares matrices in transposed form
  1418. lhs = van.T
  1419. rhs = y.T
  1420. if w is not None:
  1421. w = np.asarray(w) + 0.0
  1422. if w.ndim != 1:
  1423. raise TypeError("expected 1D vector for w")
  1424. if len(x) != len(w):
  1425. raise TypeError("expected x and w to have same length")
  1426. # apply weights. Don't use inplace operations as they
  1427. # can cause problems with NA.
  1428. lhs = lhs * w
  1429. rhs = rhs * w
  1430. # set rcond
  1431. if rcond is None:
  1432. rcond = len(x)*np.finfo(x.dtype).eps
  1433. # Determine the norms of the design matrix columns.
  1434. if issubclass(lhs.dtype.type, np.complexfloating):
  1435. scl = np.sqrt((np.square(lhs.real) + np.square(lhs.imag)).sum(1))
  1436. else:
  1437. scl = np.sqrt(np.square(lhs).sum(1))
  1438. scl[scl == 0] = 1
  1439. # Solve the least squares problem.
  1440. c, resids, rank, s = la.lstsq(lhs.T/scl, rhs.T, rcond)
  1441. c = (c.T/scl).T
  1442. # Expand c to include non-fitted coefficients which are set to zero
  1443. if deg.ndim > 0:
  1444. if c.ndim == 2:
  1445. cc = np.zeros((lmax + 1, c.shape[1]), dtype=c.dtype)
  1446. else:
  1447. cc = np.zeros(lmax + 1, dtype=c.dtype)
  1448. cc[deg] = c
  1449. c = cc
  1450. # warn on rank reduction
  1451. if rank != order and not full:
  1452. msg = "The fit may be poorly conditioned"
  1453. warnings.warn(msg, pu.RankWarning)
  1454. if full:
  1455. return c, [resids, rank, s, rcond]
  1456. else:
  1457. return c
  1458. def chebcompanion(c):
  1459. """Return the scaled companion matrix of c.
  1460. The basis polynomials are scaled so that the companion matrix is
  1461. symmetric when `c` is a Chebyshev basis polynomial. This provides
  1462. better eigenvalue estimates than the unscaled case and for basis
  1463. polynomials the eigenvalues are guaranteed to be real if
  1464. `numpy.linalg.eigvalsh` is used to obtain them.
  1465. Parameters
  1466. ----------
  1467. c : array_like
  1468. 1-D array of Chebyshev series coefficients ordered from low to high
  1469. degree.
  1470. Returns
  1471. -------
  1472. mat : ndarray
  1473. Scaled companion matrix of dimensions (deg, deg).
  1474. Notes
  1475. -----
  1476. .. versionadded::1.7.0
  1477. """
  1478. # c is a trimmed copy
  1479. [c] = pu.as_series([c])
  1480. if len(c) < 2:
  1481. raise ValueError('Series must have maximum degree of at least 1.')
  1482. if len(c) == 2:
  1483. return np.array([[-c[0]/c[1]]])
  1484. n = len(c) - 1
  1485. mat = np.zeros((n, n), dtype=c.dtype)
  1486. scl = np.array([1.] + [np.sqrt(.5)]*(n-1))
  1487. top = mat.reshape(-1)[1::n+1]
  1488. bot = mat.reshape(-1)[n::n+1]
  1489. top[0] = np.sqrt(.5)
  1490. top[1:] = 1/2
  1491. bot[...] = top
  1492. mat[:, -1] -= (c[:-1]/c[-1])*(scl/scl[-1])*.5
  1493. return mat
  1494. def chebroots(c):
  1495. """
  1496. Compute the roots of a Chebyshev series.
  1497. Return the roots (a.k.a. "zeros") of the polynomial
  1498. .. math:: p(x) = \\sum_i c[i] * T_i(x).
  1499. Parameters
  1500. ----------
  1501. c : 1-D array_like
  1502. 1-D array of coefficients.
  1503. Returns
  1504. -------
  1505. out : ndarray
  1506. Array of the roots of the series. If all the roots are real,
  1507. then `out` is also real, otherwise it is complex.
  1508. See Also
  1509. --------
  1510. polyroots, legroots, lagroots, hermroots, hermeroots
  1511. Notes
  1512. -----
  1513. The root estimates are obtained as the eigenvalues of the companion
  1514. matrix, Roots far from the origin of the complex plane may have large
  1515. errors due to the numerical instability of the series for such
  1516. values. Roots with multiplicity greater than 1 will also show larger
  1517. errors as the value of the series near such points is relatively
  1518. insensitive to errors in the roots. Isolated roots near the origin can
  1519. be improved by a few iterations of Newton's method.
  1520. The Chebyshev series basis polynomials aren't powers of `x` so the
  1521. results of this function may seem unintuitive.
  1522. Examples
  1523. --------
  1524. >>> import numpy.polynomial.chebyshev as cheb
  1525. >>> cheb.chebroots((-1, 1,-1, 1)) # T3 - T2 + T1 - T0 has real roots
  1526. array([ -5.00000000e-01, 2.60860684e-17, 1.00000000e+00])
  1527. """
  1528. # c is a trimmed copy
  1529. [c] = pu.as_series([c])
  1530. if len(c) < 2:
  1531. return np.array([], dtype=c.dtype)
  1532. if len(c) == 2:
  1533. return np.array([-c[0]/c[1]])
  1534. m = chebcompanion(c)
  1535. r = la.eigvals(m)
  1536. r.sort()
  1537. return r
  1538. def chebgauss(deg):
  1539. """
  1540. Gauss-Chebyshev quadrature.
  1541. Computes the sample points and weights for Gauss-Chebyshev quadrature.
  1542. These sample points and weights will correctly integrate polynomials of
  1543. degree :math:`2*deg - 1` or less over the interval :math:`[-1, 1]` with
  1544. the weight function :math:`f(x) = 1/\sqrt{1 - x^2}`.
  1545. Parameters
  1546. ----------
  1547. deg : int
  1548. Number of sample points and weights. It must be >= 1.
  1549. Returns
  1550. -------
  1551. x : ndarray
  1552. 1-D ndarray containing the sample points.
  1553. y : ndarray
  1554. 1-D ndarray containing the weights.
  1555. Notes
  1556. -----
  1557. .. versionadded:: 1.7.0
  1558. The results have only been tested up to degree 100, higher degrees may
  1559. be problematic. For Gauss-Chebyshev there are closed form solutions for
  1560. the sample points and weights. If n = `deg`, then
  1561. .. math:: x_i = \cos(\pi (2 i - 1) / (2 n))
  1562. .. math:: w_i = \pi / n
  1563. """
  1564. ideg = int(deg)
  1565. if ideg != deg or ideg < 1:
  1566. raise ValueError("deg must be a non-negative integer")
  1567. x = np.cos(np.pi * np.arange(1, 2*ideg, 2) / (2.0*ideg))
  1568. w = np.ones(ideg)*(np.pi/ideg)
  1569. return x, w
  1570. def chebweight(x):
  1571. """
  1572. The weight function of the Chebyshev polynomials.
  1573. The weight function is :math:`1/\sqrt{1 - x^2}` and the interval of
  1574. integration is :math:`[-1, 1]`. The Chebyshev polynomials are
  1575. orthogonal, but not normalized, with respect to this weight function.
  1576. Parameters
  1577. ----------
  1578. x : array_like
  1579. Values at which the weight function will be computed.
  1580. Returns
  1581. -------
  1582. w : ndarray
  1583. The weight function at `x`.
  1584. Notes
  1585. -----
  1586. .. versionadded:: 1.7.0
  1587. """
  1588. w = 1./(np.sqrt(1. + x) * np.sqrt(1. - x))
  1589. return w
  1590. def chebpts1(npts):
  1591. """
  1592. Chebyshev points of the first kind.
  1593. The Chebyshev points of the first kind are the points ``cos(x)``,
  1594. where ``x = [pi*(k + .5)/npts for k in range(npts)]``.
  1595. Parameters
  1596. ----------
  1597. npts : int
  1598. Number of sample points desired.
  1599. Returns
  1600. -------
  1601. pts : ndarray
  1602. The Chebyshev points of the first kind.
  1603. See Also
  1604. --------
  1605. chebpts2
  1606. Notes
  1607. -----
  1608. .. versionadded:: 1.5.0
  1609. """
  1610. _npts = int(npts)
  1611. if _npts != npts:
  1612. raise ValueError("npts must be integer")
  1613. if _npts < 1:
  1614. raise ValueError("npts must be >= 1")
  1615. x = np.linspace(-np.pi, 0, _npts, endpoint=False) + np.pi/(2*_npts)
  1616. return np.cos(x)
  1617. def chebpts2(npts):
  1618. """
  1619. Chebyshev points of the second kind.
  1620. The Chebyshev points of the second kind are the points ``cos(x)``,
  1621. where ``x = [pi*k/(npts - 1) for k in range(npts)]``.
  1622. Parameters
  1623. ----------
  1624. npts : int
  1625. Number of sample points desired.
  1626. Returns
  1627. -------
  1628. pts : ndarray
  1629. The Chebyshev points of the second kind.
  1630. Notes
  1631. -----
  1632. .. versionadded:: 1.5.0
  1633. """
  1634. _npts = int(npts)
  1635. if _npts != npts:
  1636. raise ValueError("npts must be integer")
  1637. if _npts < 2:
  1638. raise ValueError("npts must be >= 2")
  1639. x = np.linspace(-np.pi, 0, _npts)
  1640. return np.cos(x)
  1641. #
  1642. # Chebyshev series class
  1643. #
  1644. class Chebyshev(ABCPolyBase):
  1645. """A Chebyshev series class.
  1646. The Chebyshev class provides the standard Python numerical methods
  1647. '+', '-', '*', '//', '%', 'divmod', '**', and '()' as well as the
  1648. methods listed below.
  1649. Parameters
  1650. ----------
  1651. coef : array_like
  1652. Chebyshev coefficients in order of increasing degree, i.e.,
  1653. ``(1, 2, 3)`` gives ``1*T_0(x) + 2*T_1(x) + 3*T_2(x)``.
  1654. domain : (2,) array_like, optional
  1655. Domain to use. The interval ``[domain[0], domain[1]]`` is mapped
  1656. to the interval ``[window[0], window[1]]`` by shifting and scaling.
  1657. The default value is [-1, 1].
  1658. window : (2,) array_like, optional
  1659. Window, see `domain` for its use. The default value is [-1, 1].
  1660. .. versionadded:: 1.6.0
  1661. """
  1662. # Virtual Functions
  1663. _add = staticmethod(chebadd)
  1664. _sub = staticmethod(chebsub)
  1665. _mul = staticmethod(chebmul)
  1666. _div = staticmethod(chebdiv)
  1667. _pow = staticmethod(chebpow)
  1668. _val = staticmethod(chebval)
  1669. _int = staticmethod(chebint)
  1670. _der = staticmethod(chebder)
  1671. _fit = staticmethod(chebfit)
  1672. _line = staticmethod(chebline)
  1673. _roots = staticmethod(chebroots)
  1674. _fromroots = staticmethod(chebfromroots)
  1675. # Virtual properties
  1676. nickname = 'cheb'
  1677. domain = np.array(chebdomain)
  1678. window = np.array(chebdomain)