chebyshev.py 60 KB

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