chebyshev.py 61 KB

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