chebyshev.py 61 KB

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