"""Module contain algorithms for low-rank models. There is almost no typechecking of any kind here, just focus on speed """ import math from scipy.linalg import svd,inv from scipy import dot,empty,eye,newaxis,zeros,sqrt,diag,\ apply_along_axis,mean,ones,randn,empty_like,outer,r_,c_,\ rand,sum,cumsum,matrix, expand_dims,minimum,where,arange has_sym=True try: from symeig import symeig except: has_sym = False def pca(a, aopt,scale='scores',mode='normal',center_axis=0): """ Principal Component Analysis. Performs PCA on given matrix and returns results in a dictionary. :Parameters: a : array Data measurement matrix, (samples x variables) aopt : int Number of components to use, aopt<=min(samples, variables) :Returns: results : dict keys -- values, T -- scores, P -- loadings, E -- residuals, lev --leverages, ssq -- sum of squares, expvar -- cumulative explained variance, aopt -- number of components used :OtherParam eters: mode : str Amount of info retained, ('fast', 'normal', 'detailed') center_axis : int Center along given axis. If neg.: no centering (-inf,..., matrix modes) :SeeAlso: - pcr : other blm - pls : other blm - lpls : other blm Notes ----- Uses kernel speed-up if m>>n or m<>> import scipy,engines >>> a=scipy.asarray([[1,2,3],[2,4,5]]) >>> dat=engines.pca(a, 2) >>> dat['expvarx'] array([0.,99.8561562, 100.]) """ m, n = a.shape assert(aopt<=min(m,n)) if center_axis>=0: a = a - expand_dims(a.mean(center_axis), center_axis) if m>(n+100) or n>(m+100): u, s, v = esvd(a, amax=None) # fixme:amax option need to work with expl.var else: u, s, vt = svd(a, 0) v = vt.T e = s**2 tol = 1e-10 eff_rank = sum(s>s[0]*tol) aopt = minimum(aopt, eff_rank) T = u*s s = s[:aopt] T = T[:,:aopt] P = v[:,:aopt] if scale=='loads': T = T/s P = P*s if mode == 'fast': return {'T':T, 'P':P, 'aopt':aopt} if mode=='detailed': E = empty((aopt, m, n)) ssq = [] lev = [] for ai in range(aopt): E[ai,:,:] = a - dot(T[:,:ai+1], P[:,:ai+1].T) ssq.append([(E[ai,:,:]**2).mean(0), (E[ai,:,:]**2).mean(1)]) if scale=='loads': lev.append([((s*T)**2).sum(1), (P**2).sum(1)]) else: lev.append([(T**2).sum(1), ((s*P)**2).sum(1)]) else: # residuals E = a - dot(T, P.T) SEP = E**2 ssq = [SEP.sum(0), SEP.sum(1)] # leverages if scale=='loads': lev = [(1./m)+(T**2).sum(1), (1./n)+((P/s)**2).sum(1)] else: lev = [(1./m)+((T/s)**2).sum(1), (1./n)+(P**2).sum(1)] # variances expvarx = r_[0, 100*e.cumsum()/e.sum()][:aopt+1] return {'T':T, 'P':P, 'E':E, 'expvarx':expvarx, 'levx':lev, 'ssqx':ssq, 'aopt':aopt, 'eigvals': e[:aopt,newaxis]} def pcr(a, b, aopt, scale='scores',mode='normal',center_axis=0): """ Principal Component Regression. Performs PCR on given matrix and returns results in a dictionary. :Parameters: a : array Data measurement matrix, (samples x variables) b : array Data response matrix, (samples x responses) aopt : int Number of components to use, aopt<=min(samples, variables) :Returns: results : dict keys -- values, T -- scores, P -- loadings, E -- residuals, levx -- leverages, ssqx -- sum of squares, expvarx -- cumulative explained variance, aopt -- number of components used :OtherParameters: mode : str Amount of info retained, ('fast', 'normal', 'detailed') center_axis : int Center along given axis. If neg.: no centering (-inf,..., matrix modes) :SeeAlso: - pca : other blm - pls : other blm - lpls : other blm Notes ----- Uses kernel speed-up if m>>n or m<>> import scipy,engines >>> a=scipy.asarray([[1,2,3],[2,4,5]]) >>> b=scipy.asarray([[1,1],[2,3]]) >>> dat=engines.pcr(a, 2) >>> dat['expvarx'] array([0.,99.8561562, 100.]) """ k, l = m_shape(b) if center_axis>=0: b = b - expand_dims(b.mean(center_axis), center_axis) dat = pca(a, aopt=aopt, scale=scale, mode=mode, center_axis=center_axis) T = dat['T'] weights = apply_along_axis(vnorm, 0, T)**2 if scale=='loads': Q = dot(b.T, T*weights) else: Q = dot(b.T, T/weights) if mode=='fast': dat.update({'Q':Q}) return dat if mode=='detailed': F = empty((aopt, k, l)) for i in range(aopt): F[i,:,:] = b - dot(T[:,:i+1], Q[:,:i+1].T) else: F = b - dot(T, Q.T) expvary = r_[0, 100*((T**2).sum(0)*(Q**2).sum(0)/(b**2).sum()).cumsum()[:aopt]] #fixme: Y-var leverages dat.update({'Q':Q, 'F':F, 'expvary':expvary}) return dat def pls(a, b, aopt=2, scale='scores', mode='normal', center_axis=-1, ab=None): """Partial Least Squares Regression. Performs PLS on given matrix and returns results in a dictionary. :Parameters: a : array Data measurement matrix, (samples x variables) b : array Data response matrix, (samples x responses) aopt : int Number of components to use, aopt<=min(samples, variables) :Returns: results : dict keys -- values, T -- scores, P -- loadings, E -- residuals, levx -- leverages, ssqx -- sum of squares, expvarx -- cumulative explained variance of descriptors, expvary -- cumulative explained variance of responses, aopt -- number of components used :OtherParameters: mode : str Amount of info retained, ('fast', 'normal', 'detailed') center_axis : int Center along given axis. If neg.: no centering (-inf,..., matrix modes) :SeeAlso: - pca : other blm - pcr : other blm - lpls : other blm Notes ----- Uses kernel speed-up if m>>n or m<>> import scipy,engines >>> a=scipy.asarray([[1,2,3],[2,4,5]]) >>> b=scipy.asarray([[1,1],[2,3]]) >>> dat=engines.pls(a, b, 2) >>> dat['expvarx'] array([0.,99.8561562, 100.]) """ m, n = m_shape(a) if ab!=None: mm, l = m_shape(ab) assert(m==mm) else: k, l = m_shape(b) if center_axis>=0: a = a - expand_dims(a.mean(center_axis), center_axis) b = b - expand_dims(b.mean(center_axis), center_axis) W = empty((n, aopt)) P = empty((n, aopt)) R = empty((n, aopt)) Q = empty((l, aopt)) T = empty((m, aopt)) B = empty((aopt, n, l)) if ab==None: ab = dot(a.T, b) for i in range(aopt): if ab.shape[1]==1: #pls 1 w = ab.reshape(n, l) w = w/vnorm(w) elif n0: for j in range(0, i, 1): r = r - dot(P[:,j].T, w)*R[:,j][:,newaxis] t = dot(a, r) tt = vnorm(t)**2 p = dot(a.T, t)/tt q = dot(r.T, ab).T/tt ab = ab - dot(p, q.T)*tt T[:,i] = t.ravel() W[:,i] = w.ravel() if mode=='fast' and i==aopt-1: if scale=='loads': tnorm = apply_along_axis(vnorm, 0, T) T = T/tnorm W = W*tnorm return {'T':T, 'W':W} P[:,i] = p.ravel() R[:,i] = r.ravel() Q[:,i] = q.ravel() B[i] = dot(R[:,:i+1], Q[:,:i+1].T) if mode=='detailed': E = empty((aopt, m, n)) F = empty((aopt, k, l)) for i in range(1, aopt+1, 1): E[i-1] = a - dot(T[:,:i], P[:,:i].T) F[i-1] = b - dot(T[:,:i], Q[:,:i].T) else: E = a - dot(T[:,:aopt], P[:,:aopt].T) F = b - dot(T[:,:aopt], Q[:,:aopt].T) if scale=='loads': tnorm = apply_along_axis(vnorm, 0, T) T = T/tnorm W = W*tnorm Q = Q*tnorm P = P*tnorm return {'B':B, 'Q':Q, 'P':P, 'T':T, 'W':W, 'R':R, 'E':E, 'F':F} def w_simpls(aat, b, aopt): """ Simpls for wide matrices. Fast pls for crossval, used in calc rmsep for wide X There is no P or W. T is normalised """ bb = b.copy() m, m = aat.shape U = empty((m, aopt)) # W T = empty((m, aopt)) H = empty((m, aopt)) # R PROJ = empty((m, aopt)) # P? for i in range(aopt): q, s, vh = svd(dot(dot(b.T, aat), b), full_matrices=0) u = dot(b, q[:,:1]) #y-factor scores U[:,i] = u.ravel() t = dot(aat, u) t = t/vnorm(t) T[:,i] = t.ravel() h = dot(aat, t) #score-weights H[:,i] = h.ravel() PROJ[:,:i+1] = dot(T[:,:i+1], inv(dot(T[:,:i+1].T, H[:,:i+1])) ) if iY :input: X : data matrix (m, n) Y : data matrix (m, l) Z : data matrix (n, o) :output: T : X-scores W : X-weights/Z-weights P : X-loadings Q : Y-loadings U : X-Y relation L : Z-scores K : Z-loads B : Regression coefficients X->Y b0: Regression coefficient intercept evx : X-explained variance evy : Y-explained variance evz : Z-explained variance :Notes: """ if mean_ctr!=None: xctr, yctr, zctr = mean_ctr X, mnX = center(X, xctr) Y, mnY = center(Y, xctr) Z, mnZ = center(Z, zctr) varX = pow(X, 2).sum() varY = pow(Y, 2).sum() varZ = pow(Z, 2).sum() m, n = X.shape k, l = Y.shape u, o = Z.shape # initialize U = empty((k, a_max)) Q = empty((l, a_max)) T = empty((m, a_max)) W = empty((n, a_max)) P = empty((n, a_max)) K = empty((o, a_max)) L = empty((u, a_max)) B = empty((a_max, n, l)) b0 = empty((a_max, m, l)) var_x = empty((a_max,)) var_y = empty((a_max,)) var_z = empty((a_max,)) for a in range(a_max): if verbose: print "\n Working on comp. %s" %a u = Y[:,:1] diff = 1 MAX_ITER = 200 lim = 1e-16 niter = 0 while (diff>lim and niterY :input: X : data matrix (m, n) Y : data matrix (m, l) :output: T : X-scores W : X-weights/Z-weights P : X-loadings Q : Y-loadings U : X-Y relation B : Regression coefficients X->Y b0: Regression coefficient intercept evx : X-explained variance evy : Y-explained variance evz : Z-explained variance :Notes: """ if ax_center>=0: mn_x = expand_dims(X.mean(ax_center), ax_center) mn_y = expand_dims(Y.mean(ax_center), ax_center) X = X - mn_x Y = Y - mn_y varX = pow(X, 2).sum() varY = pow(Y, 2).sum() m, n = X.shape k, l = Y.shape # initialize U = empty((k, a_max)) Q = empty((l, a_max)) T = empty((m, a_max)) W = empty((n, a_max)) P = empty((n, a_max)) B = empty((a_max, n, l)) b0 = empty((a_max, m, l)) var_x = empty((a_max,)) var_y = empty((a_max,)) t1 = X[:,:1] for a in range(a_max): if verbose: print "\n Working on comp. %s" %a u = Y[:,:1] diff = 1 MAX_ITER = 100 lim = 1e-16 niter = 0 while (diff>lim and niter=n: kernel = dot(data.T, data) if has_sym: if amax==None: amax = n pcrange = None else: pcrange = [n-amax, n] s, v = symeig(kernel, range=pcrange, overwrite=True) s = s[::-1].real v = v[:,::-1].real else: u, s, vt = svd(kernel) v = vt.T s = sqrt(s) u = dot(data, v)/s else: kernel = dot(data, data.T) if has_sym: if amax==None: amax = m pcrange = None else: pcrange = [m-amax, m] s, u = symeig(kernel, range=pcrange, overwrite=True) s = s[::-1] u = u[:,::-1] else: u, s, vt = svd(kernel) s = sqrt(s) v = dot(data.T, u)/s # some use of symeig returns the 0 imaginary part return u.real, s.real, v.real def vnorm(x): # assume column arrays (or vectors) return math.sqrt(dot(x.T, x)) def center(a, axis): # 0 = col center, 1 = row center, 2 = double center # -1 = nothing if axis==-1: mn = zeros((a.shape[1],)) elif axis==0: mn = a.mean(0) elif axis==1: mn = a.mean(1)[:,newaxis] elif axis==2: mn = a.mean(0) + a.mean(1)[:,newaxis] - a.mean() else: raise IOError("input error: axis must be in [-1,0,1,2]") return a - mn, mn def scale(a, axis): if axis==-1: sc = zeros((a.shape[1],)) elif axis==0: sc = a.std(0) elif axis==1: sc = a.std(1)[:,newaxis] else: raise IOError("input error: axis must be in [-1,0,1]") return a - sc, sc