summaryrefslogtreecommitdiff
path: root/numpy/oldnumeric/mlab.py
diff options
context:
space:
mode:
authorTravis Oliphant <oliphant@enthought.com>2006-08-13 23:08:42 +0000
committerTravis Oliphant <oliphant@enthought.com>2006-08-13 23:08:42 +0000
commit9ad90f5e7264e1f1fd9ee38668c5b1e19982e9da (patch)
tree4a26a8f7fefe5fd97dc5c30bc646d9eae6b9dee8 /numpy/oldnumeric/mlab.py
parent7aa30aee6bae9a2c28539b06331225efbd3f5370 (diff)
downloadnumpy-9ad90f5e7264e1f1fd9ee38668c5b1e19982e9da.tar.gz
Restore numpy.oldnumeric.mlab.cov to MLab.cov behavior
Diffstat (limited to 'numpy/oldnumeric/mlab.py')
-rw-r--r--numpy/oldnumeric/mlab.py33
1 files changed, 29 insertions, 4 deletions
diff --git a/numpy/oldnumeric/mlab.py b/numpy/oldnumeric/mlab.py
index 3d39ef0a9..216ef1391 100644
--- a/numpy/oldnumeric/mlab.py
+++ b/numpy/oldnumeric/mlab.py
@@ -8,7 +8,7 @@ from numpy import tril, trapz as _Ntrapz, hanning, rot90, triu, diff, \
angle, roots, ptp as _Nptp, kaiser, cumprod as _Ncumprod, \
diag, msort, prod as _Nprod, std as _Nstd, hamming, flipud, \
amax as _Nmax, amin as _Nmin, blackman, bartlett, corrcoef as _Ncorrcoef,\
- cov as _Ncov, squeeze, sinc, median, fliplr, mean as _Nmean
+ cov as _Ncov, squeeze, sinc, median, fliplr, mean as _Nmean, transpose
from numpy.linalg import eig, svd
from numpy.random import rand, randn
@@ -59,11 +59,36 @@ def std(x, axis=0):
def mean(x, axis=0):
return _Nmean(x, axis)
+# This is exactly the same cov function as in MLab
def cov(m, y=None, rowvar=0, bias=0):
- return _Ncov(m, y, rowvar, bias)
-
+ if y is None:
+ y = m
+ else:
+ y = y
+ if rowvar:
+ m = transpose(m)
+ y = transpose(y)
+ if (m.shape[0] == 1):
+ m = transpose(m)
+ if (y.shape[0] == 1):
+ y = transpose(y)
+ N = m.shape[0]
+ if (y.shape[0] != N):
+ raise ValueError, "x and y must have the same number "\
+ "of observations"
+ m = m - _Nmean(m,axis=0)
+ y = y - _Nmean(y,axis=0)
+ if bias:
+ fact = N*1.0
+ else:
+ fact = N-1.0
+ return squeeze(dot(transpose(m), conjugate(y)) / fact)
+
+from numpy import sqrt, multiply
def corrcoef(x, y=None):
- return _Ncorrcoef(x,y,0,0)
+ c = cov(x,y)
+ d = diag(c)
+ return c/sqrt(multiply.outer(d,d))
from compat import *
from functions import *