1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464
|
"""Tests for chebyshev module.
"""
from __future__ import division
import numpy as np
import numpy.polynomial.chebyshev as ch
from numpy.testing import *
from exceptions import TypeError, ValueError
def trim(x) :
return ch.chebtrim(x, tol=1e-6)
T0 = [ 1]
T1 = [ 0, 1]
T2 = [-1, 0, 2]
T3 = [ 0, -3, 0, 4]
T4 = [ 1, 0, -8, 0, 8]
T5 = [ 0, 5, 0, -20, 0, 16]
T6 = [-1, 0, 18, 0, -48, 0, 32]
T7 = [ 0, -7, 0, 56, 0, -112, 0, 64]
T8 = [ 1, 0, -32, 0, 160, 0, -256, 0, 128]
T9 = [ 0, 9, 0, -120, 0, 432, 0, -576, 0, 256]
Tlist = [T0, T1, T2, T3, T4, T5, T6, T7, T8, T9]
class TestPrivate(TestCase) :
def test__cseries_to_zseries(self) :
for i in range(5) :
inp = np.array([2] + [1]*i, np.double)
tgt = np.array([.5]*i + [2] + [.5]*i, np.double)
res = ch._cseries_to_zseries(inp)
assert_equal(res, tgt)
def test__zseries_to_cseries(self) :
for i in range(5) :
inp = np.array([.5]*i + [2] + [.5]*i, np.double)
tgt = np.array([2] + [1]*i, np.double)
res = ch._zseries_to_cseries(inp)
assert_equal(res, tgt)
class TestConstants(TestCase) :
def test_chebdomain(self) :
assert_equal(ch.chebdomain, [-1, 1])
def test_chebzero(self) :
assert_equal(ch.chebzero, [0])
def test_chebone(self) :
assert_equal(ch.chebone, [1])
def test_chebx(self) :
assert_equal(ch.chebx, [0, 1])
class TestArithmetic(TestCase) :
def test_chebadd(self) :
for i in range(5) :
for j in range(5) :
msg = "At i=%d, j=%d" % (i,j)
tgt = np.zeros(max(i,j) + 1)
tgt[i] += 1
tgt[j] += 1
res = ch.chebadd([0]*i + [1], [0]*j + [1])
assert_equal(trim(res), trim(tgt), err_msg=msg)
def test_chebsub(self) :
for i in range(5) :
for j in range(5) :
msg = "At i=%d, j=%d" % (i,j)
tgt = np.zeros(max(i,j) + 1)
tgt[i] += 1
tgt[j] -= 1
res = ch.chebsub([0]*i + [1], [0]*j + [1])
assert_equal(trim(res), trim(tgt), err_msg=msg)
def test_chebmul(self) :
for i in range(5) :
for j in range(5) :
msg = "At i=%d, j=%d" % (i,j)
tgt = np.zeros(i + j + 1)
tgt[i + j] += .5
tgt[abs(i - j)] += .5
res = ch.chebmul([0]*i + [1], [0]*j + [1])
assert_equal(trim(res), trim(tgt), err_msg=msg)
def test_chebdiv(self) :
for i in range(5) :
for j in range(5) :
msg = "At i=%d, j=%d" % (i,j)
ci = [0]*i + [1]
cj = [0]*j + [1]
tgt = ch.chebadd(ci, cj)
quo, rem = ch.chebdiv(tgt, ci)
res = ch.chebadd(ch.chebmul(quo, ci), rem)
assert_equal(trim(res), trim(tgt), err_msg=msg)
def test_chebval(self) :
def f(x) :
return x*(x**2 - 1)
#check empty input
assert_equal(ch.chebval([], [1]).size, 0)
#check normal input)
for i in range(5) :
tgt = 1
res = ch.chebval(1, [0]*i + [1])
assert_almost_equal(res, tgt)
tgt = (-1)**i
res = ch.chebval(-1, [0]*i + [1])
assert_almost_equal(res, tgt)
zeros = np.cos(np.linspace(-np.pi, 0, 2*i + 1)[1::2])
tgt = 0
res = ch.chebval(zeros, [0]*i + [1])
assert_almost_equal(res, tgt)
x = np.linspace(-1,1)
tgt = f(x)
res = ch.chebval(x, [0, -.25, 0, .25])
assert_almost_equal(res, tgt)
#check that shape is preserved
for i in range(3) :
dims = [2]*i
x = np.zeros(dims)
assert_equal(ch.chebval(x, [1]).shape, dims)
assert_equal(ch.chebval(x, [1,0]).shape, dims)
assert_equal(ch.chebval(x, [1,0,0]).shape, dims)
class TestCalculus(TestCase) :
def test_chebint(self) :
# check exceptions
assert_raises(ValueError, ch.chebint, [0], -1)
assert_raises(ValueError, ch.chebint, [0], 1, [0,0])
assert_raises(ValueError, ch.chebint, [0], 1, lbnd=[0,0])
assert_raises(ValueError, ch.chebint, [0], 1, scl=[0,0])
# check single integration with integration constant
for i in range(5) :
scl = i + 1
pol = [0]*i + [1]
tgt = [i] + [0]*i + [1/scl]
chebpol = ch.poly2cheb(pol)
chebint = ch.chebint(chebpol, m=1, k=[i])
res = ch.cheb2poly(chebint)
assert_almost_equal(trim(res), trim(tgt))
# check single integration with integration constant and lbnd
for i in range(5) :
scl = i + 1
pol = [0]*i + [1]
chebpol = ch.poly2cheb(pol)
chebint = ch.chebint(chebpol, m=1, k=[i], lbnd=-1)
assert_almost_equal(ch.chebval(-1, chebint), i)
# check single integration with integration constant and scaling
for i in range(5) :
scl = i + 1
pol = [0]*i + [1]
tgt = [i] + [0]*i + [2/scl]
chebpol = ch.poly2cheb(pol)
chebint = ch.chebint(chebpol, m=1, k=[i], scl=2)
res = ch.cheb2poly(chebint)
assert_almost_equal(trim(res), trim(tgt))
# check multiple integrations with default k
for i in range(5) :
for j in range(2,5) :
pol = [0]*i + [1]
tgt = pol[:]
for k in range(j) :
tgt = ch.chebint(tgt, m=1)
res = ch.chebint(pol, m=j)
assert_almost_equal(trim(res), trim(tgt))
# check multiple integrations with defined k
for i in range(5) :
for j in range(2,5) :
pol = [0]*i + [1]
tgt = pol[:]
for k in range(j) :
tgt = ch.chebint(tgt, m=1, k=[k])
res = ch.chebint(pol, m=j, k=range(j))
assert_almost_equal(trim(res), trim(tgt))
# check multiple integrations with lbnd
for i in range(5) :
for j in range(2,5) :
pol = [0]*i + [1]
tgt = pol[:]
for k in range(j) :
tgt = ch.chebint(tgt, m=1, k=[k], lbnd=-1)
res = ch.chebint(pol, m=j, k=range(j), lbnd=-1)
assert_almost_equal(trim(res), trim(tgt))
# check multiple integrations with scaling
for i in range(5) :
for j in range(2,5) :
pol = [0]*i + [1]
tgt = pol[:]
for k in range(j) :
tgt = ch.chebint(tgt, m=1, k=[k], scl=2)
res = ch.chebint(pol, m=j, k=range(j), scl=2)
assert_almost_equal(trim(res), trim(tgt))
def test_chebder(self) :
# check exceptions
assert_raises(ValueError, ch.chebder, [0], -1)
# check that zeroth deriviative does nothing
for i in range(5) :
tgt = [1] + [0]*i
res = ch.chebder(tgt, m=0)
assert_equal(trim(res), trim(tgt))
# check that derivation is the inverse of integration
for i in range(5) :
for j in range(2,5) :
tgt = [1] + [0]*i
res = ch.chebder(ch.chebint(tgt, m=j), m=j)
assert_almost_equal(trim(res), trim(tgt))
# check derivation with scaling
for i in range(5) :
for j in range(2,5) :
tgt = [1] + [0]*i
res = ch.chebder(ch.chebint(tgt, m=j, scl=2), m=j, scl=.5)
assert_almost_equal(trim(res), trim(tgt))
class TestMisc(TestCase) :
def test_chebfromroots(self) :
res = ch.chebfromroots([])
assert_almost_equal(trim(res), [1])
for i in range(1,5) :
roots = np.cos(np.linspace(-np.pi, 0, 2*i + 1)[1::2])
tgt = [0]*i + [1]
res = ch.chebfromroots(roots)*2**(i-1)
assert_almost_equal(trim(res),trim(tgt))
def test_chebroots(self) :
assert_almost_equal(ch.chebroots([1]), [])
assert_almost_equal(ch.chebroots([1, 2]), [-.5])
for i in range(2,5) :
tgt = np.linspace(-1, 1, i)
res = ch.chebroots(ch.chebfromroots(tgt))
assert_almost_equal(trim(res), trim(tgt))
def test_chebvander(self) :
# check for 1d x
x = np.arange(3)
v = ch.chebvander(x, 3)
assert_(v.shape == (3,4))
for i in range(4) :
coef = [0]*i + [1]
assert_almost_equal(v[...,i], ch.chebval(x, coef))
# check for 2d x
x = np.array([[1,2],[3,4],[5,6]])
v = ch.chebvander(x, 3)
assert_(v.shape == (3,2,4))
for i in range(4) :
coef = [0]*i + [1]
assert_almost_equal(v[...,i], ch.chebval(x, coef))
def test_chebfit(self) :
def f(x) :
return x*(x - 1)*(x - 2)
# Test exceptions
assert_raises(ValueError, ch.chebfit, [1], [1], -1)
assert_raises(TypeError, ch.chebfit, [[1]], [1], 0)
assert_raises(TypeError, ch.chebfit, [], [1], 0)
assert_raises(TypeError, ch.chebfit, [1], [[[1]]], 0)
assert_raises(TypeError, ch.chebfit, [1, 2], [1], 0)
assert_raises(TypeError, ch.chebfit, [1], [1, 2], 0)
# Test fit
x = np.linspace(0,2)
y = f(x)
coef = ch.chebfit(x, y, 3)
assert_equal(len(coef), 4)
assert_almost_equal(ch.chebval(x, coef), y)
coef = ch.chebfit(x, y, 4)
assert_equal(len(coef), 5)
assert_almost_equal(ch.chebval(x, coef), y)
coef2d = ch.chebfit(x, np.array([y,y]).T, 4)
assert_almost_equal(coef2d, np.array([coef,coef]).T)
def test_chebtrim(self) :
coef = [2, -1, 1, 0]
# Test exceptions
assert_raises(ValueError, ch.chebtrim, coef, -1)
# Test results
assert_equal(ch.chebtrim(coef), coef[:-1])
assert_equal(ch.chebtrim(coef, 1), coef[:-3])
assert_equal(ch.chebtrim(coef, 2), [0])
def test_chebline(self) :
assert_equal(ch.chebline(3,4), [3, 4])
def test_cheb2poly(self) :
for i in range(10) :
assert_equal(ch.cheb2poly([0]*i + [1]), Tlist[i])
def test_poly2cheb(self) :
for i in range(10) :
assert_equal(ch.poly2cheb(Tlist[i]), [0]*i + [1])
class TestChebyshevClass(TestCase) :
p1 = ch.Chebyshev([1,2,3])
p2 = ch.Chebyshev([1,2,3], [0,1])
p3 = ch.Chebyshev([1,2])
p4 = ch.Chebyshev([2,2,3])
p5 = ch.Chebyshev([3,2,3])
def test_equal(self) :
assert_(self.p1 == self.p1)
assert_(self.p2 == self.p2)
assert_(not self.p1 == self.p2)
assert_(not self.p1 == self.p3)
assert_(not self.p1 == [1,2,3])
def test_not_equal(self) :
assert_(not self.p1 != self.p1)
assert_(not self.p2 != self.p2)
assert_(self.p1 != self.p2)
assert_(self.p1 != self.p3)
assert_(self.p1 != [1,2,3])
def test_add(self) :
tgt = ch.Chebyshev([2,4,6])
assert_(self.p1 + self.p1 == tgt)
assert_(self.p1 + [1,2,3] == tgt)
assert_([1,2,3] + self.p1 == tgt)
def test_sub(self) :
tgt = ch.Chebyshev([1])
assert_(self.p4 - self.p1 == tgt)
assert_(self.p4 - [1,2,3] == tgt)
assert_([2,2,3] - self.p1 == tgt)
def test_mul(self) :
tgt = ch.Chebyshev([7.5, 10., 8., 6., 4.5])
assert_(self.p1 * self.p1 == tgt)
assert_(self.p1 * [1,2,3] == tgt)
assert_([1,2,3] * self.p1 == tgt)
def test_floordiv(self) :
tgt = ch.Chebyshev([1])
assert_(self.p4 // self.p1 == tgt)
assert_(self.p4 // [1,2,3] == tgt)
assert_([2,2,3] // self.p1 == tgt)
def test_mod(self) :
tgt = ch.Chebyshev([1])
assert_((self.p4 % self.p1) == tgt)
assert_((self.p4 % [1,2,3]) == tgt)
assert_(([2,2,3] % self.p1) == tgt)
def test_divmod(self) :
tquo = ch.Chebyshev([1])
trem = ch.Chebyshev([2])
quo, rem = divmod(self.p5, self.p1)
assert_(quo == tquo and rem == trem)
quo, rem = divmod(self.p5, [1,2,3])
assert_(quo == tquo and rem == trem)
quo, rem = divmod([3,2,3], self.p1)
assert_(quo == tquo and rem == trem)
def test_pow(self) :
tgt = ch.Chebyshev([1])
for i in range(5) :
res = self.p1**i
assert_(res == tgt)
tgt *= self.p1
def test_call(self) :
# domain = [-1, 1]
x = np.linspace(-1, 1)
tgt = 3*(2*x**2 - 1) + 2*x + 1
assert_almost_equal(self.p1(x), tgt)
# domain = [0, 1]
x = np.linspace(0, 1)
xx = 2*x - 1
assert_almost_equal(self.p2(x), self.p1(xx))
def test_convert(self) :
x = np.linspace(-1,1)
p = self.p1.convert(domain=[0,1])
assert_almost_equal(p(x), self.p1(x))
def test_mapparms(self) :
parms = self.p2.mapparms()
assert_almost_equal(parms, [-1, 2])
def test_trim(self) :
coef = [1, 1e-6, 1e-12, 0]
p = ch.Chebyshev(coef)
assert_equal(p.trim().coef, coef[:3])
assert_equal(p.trim(1e-10).coef, coef[:2])
assert_equal(p.trim(1e-5).coef, coef[:1])
def test_truncate(self) :
assert_raises(ValueError, self.p1.truncate, 0)
assert_equal(len(self.p1.truncate(4)), 3)
assert_equal(len(self.p1.truncate(3)), 3)
assert_equal(len(self.p1.truncate(2)), 2)
assert_equal(len(self.p1.truncate(1)), 1)
def test_copy(self) :
p = self.p1.copy()
assert_(self.p1 == p)
def test_integ(self) :
p = self.p2.integ()
assert_almost_equal(p.coef, ch.chebint([1,2,3], 1, 0, scl=.5))
p = self.p2.integ(lbnd=0)
assert_almost_equal(p(0), 0)
p = self.p2.integ(1, 1)
assert_almost_equal(p.coef, ch.chebint([1,2,3], 1, 1, scl=.5))
p = self.p2.integ(2, [1, 2])
assert_almost_equal(p.coef, ch.chebint([1,2,3], 2, [1,2], scl=.5))
def test_deriv(self) :
p = self.p2.integ(2, [1, 2])
assert_almost_equal(p.deriv(1).coef, self.p2.integ(1, [1]).coef)
assert_almost_equal(p.deriv(2).coef, self.p2.coef)
def test_roots(self) :
p = ch.Chebyshev(ch.poly2cheb([0, -1, 0, 1]), [0, 1])
res = p.roots()
tgt = [0, .5, 1]
assert_almost_equal(res, tgt)
def test_fromroots(self) :
roots = [0, .5, 1]
p = ch.Chebyshev.fromroots(roots, domain=[0, 1])
res = p.coef
tgt = ch.poly2cheb([0, -1, 0, 1])
assert_almost_equal(res, tgt)
def test_fit(self) :
def f(x) :
return x*(x - 1)*(x - 2)
x = np.linspace(0,3)
y = f(x)
p = ch.Chebyshev.fit(x, y, 3)
assert_almost_equal(p(x), y)
p = ch.Chebyshev.fit(x, y, 3, None)
assert_almost_equal(p(x), y)
assert_almost_equal(p.domain, [0,3])
def test_identity(self) :
x = np.linspace(0,3)
p = ch.Chebyshev.identity()
assert_almost_equal(p(x), x)
p = ch.Chebyshev.identity([1,3])
assert_almost_equal(p(x), x)
|