Package csb :: Package numeric
[frames] | no frames]

Source Code for Package csb.numeric

  1  """ 
  2  Low level numeric / math utility functions. 
  3  """ 
  4   
  5  import sys 
  6  import math 
  7  import numpy 
  8   
  9   
 10  EXP_MIN = -308 
 11  EXP_MAX = +709 
 12       
 13  LOG_MIN = 1e-308 
 14  LOG_MAX = 1e+308 
 15   
 16  ## Euler-Mascheroni constant 
 17  EULER_MASCHERONI = 0.57721566490153286060651209008240243104215933593992 
 18   
 19   
20 -def log(x, x_min=LOG_MIN, x_max=LOG_MAX):
21 """ 22 Safe version of log, clips argument such that overflow does not occur. 23 24 @param x: input 25 @type x: numpy array or float or int 26 27 @param x_min: lower value for clipping 28 @type x_min: float 29 30 @param x_max: upper value for clipping 31 @type x_max: float 32 """ 33 34 x_min = max(x_min, LOG_MIN) 35 x_max = min(x_max, LOG_MAX) 36 37 return numpy.log(numpy.clip(x, x_min, x_max))
38
39 -def exp(x, x_min=EXP_MIN, x_max=EXP_MAX):
40 """ 41 Safe version of exp, clips argument such that overflow does not occur. 42 43 @param x: input 44 @type x: numpy array or float or int 45 46 @param x_min: lower value for clipping 47 @type x_min: float 48 49 @param x_max: upper value for clipping 50 @type x_max: float 51 """ 52 53 x_min = max(x_min, EXP_MIN) 54 x_max = min(x_max, EXP_MAX) 55 56 return numpy.exp(numpy.clip(x, x_min, x_max))
57
58 -def sign(x):
59 """ 60 Return the sign of the input. 61 """ 62 return numpy.sign(x)
63
64 -def isreal(x, tol=1e-14):
65 """ 66 Check if input array has no imaginary part. 67 68 @param x: input array 69 @type x: numpy array 70 71 @param tol: tolerance to check for equality zero 72 @type tol: float 73 """ 74 return not hasattr(x, 'real') or abs(x.imag) < tol
75
76 -def log_sum_exp(x, axis=0):
77 """ 78 Return the logarithm of the sum of exponentials. 79 80 @type x: Numpy array 81 """ 82 xmax = x.max(axis) 83 return log(exp(x - xmax).sum(axis)) + xmax
84
85 -def log_sum_exp_accumulate(x, axis=0):
86 """ 87 Return the logarithm of the accumulated sums of exponentials. 88 89 @type x: Numpy array 90 """ 91 xmax = x.max(axis) 92 return log(numpy.add.accumulate(exp(x - xmax), axis)) + xmax
93
94 -def radian2degree(x):
95 """ 96 Convert radians angles to torsion angles. 97 98 @param x: radian angle 99 @return: torsion angle of x 100 """ 101 x = x % (2 * numpy.pi) 102 numpy.putmask(x, x > numpy.pi, x - 2 * numpy.pi) 103 return x * 180. / numpy.pi
104
105 -def degree2radian(x):
106 """ 107 Convert randian angles to torsion angles. 108 109 @param x: torsion angle 110 @return: radian angle of x 111 """ 112 numpy.putmask(x, x < 0., x + 360.) 113 return x * numpy.pi / 180.
114
115 -def euler_angles(r):
116 """ 117 Calculate the euler angles from a three dimensional rotation matrix. 118 119 @param r: 3x3 Rotation matrix 120 """ 121 a = numpy.arctan2(r[2, 1], r[2, 0]) % (2 * numpy.pi) 122 b = numpy.arctan2((r[2, 0] + r[2, 1]) / (numpy.cos(a) + numpy.sin(a)), r[2, 2]) % (2 * numpy.pi) 123 c = numpy.arctan2(r[1, 2], -r[0, 2]) % (2 * numpy.pi) 124 125 return a, b, c
126
127 -def euler(a, b, c):
128 """ 129 Calculate a three dimensional rotation matrix from the euler angles. 130 131 @param a: alpha, angle between the x-axis and the line of nodes 132 @param b: beta, angle between the z axis of the different coordinate systems 133 @param c: gamma, angle between the line of nodes and the X-axis 134 """ 135 from numpy import cos, sin, array 136 137 ca, cb, cc = cos(a), cos(b), cos(c) 138 sa, sb, sc = sin(a), sin(b), sin(c) 139 140 return array([[ cc * cb * ca - sc * sa, cc * cb * sa + sc * ca, -cc * sb], 141 [-sc * cb * ca - cc * sa, -sc * cb * sa + cc * ca, sc * sb], 142 [ sb * ca, sb * sa, cb ]])
143
144 -def norm(x):
145 """ 146 Calculate the Eucledian norm of a d-dimensional vector. 147 148 @param x: vector (i.e. rank one array) 149 @return: length of vector 150 """ 151 return numpy.linalg.norm(x)
152
153 -def reverse(array, axis=0):
154 """ 155 Reverse the order of elements in an array. 156 """ 157 from numpy import take, arange 158 return take(array, arange(array.shape[axis] - 1, -1, -1), axis)
159
160 -def polar(x):
161 """ 162 Polar coordinate representation of a d-dimensional vector. 163 164 @param x: vector (i.e. rank one array) 165 @return: polar coordinates (radius and polar angles) 166 """ 167 168 (d,) = x.shape 169 phi = numpy.zeros(d) 170 for i in reversed(range(1, d)): 171 phi[i - 1] = numpy.arctan2(x[i] / numpy.cos(phi[i]), x[i - 1]) 172 173 return numpy.array([norm(x)] + phi[:-1].tolist())
174
175 -def from_polar(x):
176 """ 177 Reconstruct d-dimensional vector from polar coordinates. 178 179 @param x: vector (i.e. rank one array) 180 @return: position in d-dimensional space 181 """ 182 183 (d,) = x.shape 184 185 c = numpy.cos(x[1:]) 186 s = numpy.sin(x[1:]) 187 r = x[0] 188 189 x = numpy.zeros(d) 190 x[0] = r 191 192 for i in range(d - 1): 193 194 x[i + 1] = x[i] * s[i] 195 x[i] *= c[i] 196 197 return x
198
199 -def polar3d(x):
200 """ 201 Polar coordinate representation of a three-dimensional vector. 202 203 @param x: vector (i.e. rank one array) 204 @return: polar coordinates (radius and polar angles) 205 """ 206 207 if x.shape != (3,): 208 raise ValueError(x) 209 210 r = norm(x) 211 theta = numpy.arccos(x[2] / r) 212 phi = numpy.arctan2(x[1], x[0]) 213 214 return numpy.array([r, theta, phi])
215
216 -def from_polar3d(x):
217 """ 218 Reconstruct 3-dimensional vector from polar coordinates. 219 220 @param x: vector (i.e. rank one array) 221 @return: position in 3-dimensional space 222 """ 223 assert x.shape == (3,) 224 225 r, theta, phi = x[:] 226 s = numpy.sin(theta) 227 c = numpy.cos(theta) 228 S = numpy.sin(phi) 229 C = numpy.cos(phi) 230 231 return r * numpy.array([s * C, s * S, c])
232
233 -def dihedral_angle(a, b, c, d):
234 """ 235 Calculate the dihedral angle between 4 vectors, 236 representing 4 connected points. The angle is in range [-180, 180]. 237 238 @param a: the four points that define the dihedral angle 239 @type a: array 240 241 @return: angle in [-180, 180] 242 """ 243 244 v = b - c 245 m = numpy.cross((a - b), v) 246 m /= norm(m) 247 n = numpy.cross((d - c), v) 248 n /= norm(m) 249 250 c = numpy.dot(m, n) 251 s = numpy.dot(numpy.cross(n, m), v) / norm(v) 252 253 angle = math.degrees(math.atan2(s, c)) 254 255 if angle > 0: 256 return numpy.fmod(angle + 180, 360) - 180 257 else: 258 return numpy.fmod(angle - 180, 360) + 180
259
260 -def psi(x):
261 """ 262 Digamma function 263 """ 264 from numpy import inf, log, sum, exp 265 266 coef = [-1. / 12., 1. / 120., -1. / 252., 1. / 240., -1. / 132., 267 691. / 32760., -1. / 12.] 268 269 if x == 0.: 270 return -inf 271 elif x < 0.: 272 raise ValueError('not defined for negative values') 273 elif x < 6.: 274 return psi(x + 1) - 1. / x 275 else: 276 logx = log(x) 277 res = logx - 0.5 / x 278 res += sum([coef[i] * exp(-2 * (i + 1) * logx) for i in range(7)]) 279 return res
280
281 -def approx_psi(x):
282 from numpy import log, clip, where 283 284 if type(x) == numpy.ndarray: 285 y = 0. * x 286 y[where(x < 0.6)] = -EULER_MASCHERONI - 1. / clip(x[where(x < 0.6)], 1e-154, 1e308) 287 y[where(x >= 0.6)] = log(x[where(x >= 0.6)] - 0.5) 288 return y 289 else: 290 if x < 0.6: 291 return -EULER_MASCHERONI - 1 / clip(x, 1e-154, 1e308) 292 else: 293 return log(x - 0.5)
294
295 -def d_approx_psi(x):
296 from numpy import clip, where 297 298 if type(x) == numpy.ndarray: 299 y = 0. * x 300 y[where(x < 0.6)] = 1. / clip(x[where(x < 0.6)], 1e-154, 1e308) ** 2 301 y[where(x >= 0.6)] = 1. / (x[where(x >= 0.6)] - 0.5) 302 return y 303 else: 304 305 if x < 0.6: 306 return 1 / clip(x, 1e-154, 1e308) ** 2 307 else: 308 return 1 / (x - 0.5)
309
310 -def inv_psi(y, tol=1e-10, n_iter=100, psi=psi):
311 """ 312 Inverse digamma function 313 """ 314 from numpy import exp 315 from scipy.special import digamma 316 ## initial value 317 318 if y < -5 / 3. - EULER_MASCHERONI: 319 x = -1 / (EULER_MASCHERONI + y) 320 else: 321 x = 0.5 + exp(y) 322 323 ## Newton root finding 324 325 for dummy in range(n_iter): 326 327 z = digamma(x) - y 328 329 if abs(z) < tol: 330 break 331 332 x -= z / d_approx_psi(x) 333 334 return x
335
336 -def log_trapezoidal(log_y, x=None):
337 """ 338 Compute the logarithm of the 1D integral of x, using trepezoidal approximation. 339 Assumes x is monotonically increasing. 340 """ 341 if x is not None: 342 log_x_diff = log(x[1:] - x[:-1]) 343 else: 344 log_x_diff = 0. 345 346 log_y_add = log_sum_exp(numpy.vstack((log_y[:-1], log_y[1:])), 0) 347 348 return log_sum_exp(log_y_add + log_x_diff) - log(2)
349
350 -def log_midpoint_rule_2d(log_f, x, y):
351 x_delta = x[:-1] - x[1:] 352 y_delta = y[:-1] - y[1:] 353 354 z = numpy.array([log_f[:, 1:] , log_f[:, :-1]]) 355 356 y_hat = log_sum_exp(z.reshape((2, -1)), 0) 357 y_hat = numpy.reshape(y_hat, (len(x), len(y) - 1)) 358 y_hat += log(y_delta) - log(2) 359 360 return log_sum_exp(y_hat + log(x_delta)) - log(2.)
361
362 -def log_trapezoidal_2d(log_f, x=None, y=None):
363 """ 364 Compute the logarithm of the 1D integral of x, using trepezoidal approximation. 365 Assumes x and y is monotonically increasing. 366 """ 367 int_y = numpy.array([log_trapezoidal(log_f[i, :], y) for i in range(len(y))]) 368 369 return log_trapezoidal(int_y, x)
370
371 -def trapezoidal(x, y):
372 return 0.5 * numpy.dot(x[1:] - x[:-1], y[1:] + y[:-1])
373
374 -def trapezoidal_2d(f):
375 """ 376 Approximate the integral of f from a to b in two dimensions, 377 using trepezoidal approximation. 378 379 @param f: 2D numpy array of function values at equally spaces positions 380 @return: approximation of the definit integral 381 """ 382 383 I = f[0, 0] + f[-1, -1] + f[0, -1] + f[-1, 0] 384 I += 2 * (f[1:-1, (0, -1)].sum() + f[(0, -1), 1:-1].sum()) 385 I += 4 * f[1:-1, 1:-1].sum() 386 387 return I / 4.
388
389 -def simpson_2d(f):
390 """ 391 Approximate the integral of f from a to b in two dimensions, 392 using Composite Simpson's rule. 393 394 @param f: 2D numpy array of function values 395 @return: approximation of the definit integral 396 """ 397 398 n = int((f.shape[0] - 1) / 2) 399 i = 2 * numpy.arange(1, n + 1) - 1 400 j = 2 * numpy.arange(1, n) 401 402 I = f[0, 0] + f[-1, -1] + f[0, -1] + f[-1, 0] 403 I += 4 * (f[0, i].sum() + f[-1, i].sum() + f[0, j].sum() + f[-1, j].sum()) 404 I += 4 * (f[i, 0].sum() + f[i, -1].sum() + f[j, 0].sum() + f[j, -1].sum()) 405 I += 16 * f[i][:, i].sum() + 8 * (f[i][:, j].sum() + f[j][:, i].sum()) 406 I += 4 * f[j][:, j].sum() 407 408 return I / 9.
409
410 -def pad(x, s):
411 """ 412 Add layers of zeros around grid. 413 """ 414 s = numpy.array(s) - 1 415 y = numpy.zeros(numpy.array(x.shape) + s) 416 s /= 2 417 slices = [slice(s[i], -s[i]) for i in range(len(s))] 418 y[slices] = x 419 420 return y
421
422 -def trim(x, s):
423 """ 424 Remove additional layers. 425 """ 426 s = numpy.array(s) - 1 427 s /= 2 428 429 slices = [slice(s[i], -s[i]) for i in range(len(s))] 430 431 return x[slices]
432
433 -def zerofill(x, s):
434 435 y = numpy.zeros(s) 436 slices = [slice(-x.shape[i], None) for i in range(len(s))] 437 y[slices] = x 438 439 return y
440
441 -def convolve(x, f):
442 443 from numpy import fft, all 444 445 sx = numpy.array(x.shape) 446 sf = numpy.array(f.shape) 447 448 if not all(sx >= sf): return convolve(f, x) 449 450 y = fft.ifftn(fft.fftn(x) * fft.fftn(f, sx)).real 451 slices = [slice(sf[i] - 1, sx[i]) for i in range(len(sf))] 452 453 return y[slices]
454
455 -def correlate(x, y):
456 457 from numpy import fft 458 459 sx = numpy.array(x.shape) 460 sy = numpy.array(y.shape) 461 462 if (sx >= sy).sum(): 463 464 slices = [slice(None, sx[i] - sy[i] + 1) for i in range(len(sx))] 465 466 X = fft.fftn(x) 467 Y = fft.fftn(zerofill(y, sx)) 468 469 else: 470 471 sf = sx + sy - 1 472 slices = [slice(None, sf[i]) for i in range(len(sf))] 473 474 X = fft.fftn(x, sf) 475 Y = fft.fftn(zerofill(y, sf), sf) 476 477 return fft.ifftn(X.conjugate() * Y)[slices].real
478 479
480 -def gower_matrix(X):
481 """ 482 Gower, J.C. (1966). Some distance properties of latent root 483 and vector methods used in multivariate analysis. 484 Biometrika 53: 325-338 485 486 @param X: ensemble coordinates 487 @type X: (m,n,k) numpy.array 488 489 @return: symmetric dissimilarity matrix 490 @rtype: (n,n) numpy.array 491 """ 492 X = numpy.asarray(X) 493 494 B = sum(numpy.dot(x, x.T) for x in X) / float(len(X)) 495 b = B.mean(1) 496 bb = b.mean() 497 498 return (B - numpy.add.outer(b, b)) + bb
499