acgc.stats.sma_warton_fit
Alternate implementation of Robust SMA line fitting
This module implements the robust fitting method of Taskinen and Warton
(2013, http://dx.doi.org/10.1016/j.jtbi.2013.05.010).
The main difference from acgc.stats.bivariate_lines.sma
is the method
for computing the robust covariance. For normal (non-robust) SMA, both
functions should give the same results.
1# -*- coding: utf-8 -*- 2"""Alternate implementation of Robust SMA line fitting 3 4This module implements the robust fitting method of Taskinen and Warton 5(2013, http://dx.doi.org/10.1016/j.jtbi.2013.05.010). 6The main difference from `acgc.stats.bivariate_lines.sma` is the method 7for computing the robust covariance. For normal (non-robust) SMA, both 8functions should give the same results. 9 10""" 11import numpy as np 12import scipy.linalg as la 13import scipy.stats as stats 14 15#from numba import jit 16 17__all__ = ['sma_warton'] 18 19def sma_warton(X0,Y0,alpha=0.05,intercept=True,robust=False,wartonrobust=False): 20 '''Standard major axis, aka reduced major axis 21 22 The main advantage of this is that the best fit 23 of Y to X will be the same as the best fit of X to Y. 24 More details in Warton et al. Biology Review 2006 25 26 Parameters 27 ---------- 28 X, Y : array_like 29 Input values, Must have same length. 30 alpha : float 31 Desired confidence level for output. 32 33 Returns 34 ------- 35 Gradient : float 36 Gradient or Slope of Y vs. X 37 Intercept : floats 38 Y intercept. 39 stdgrad : float 40 Standard error of gradient estimate 41 stdint : float 42 standard error of intercept estimate 43 ci_grad : [float, float] 44 confidence interval for gradient at confidence level alpha 45 ci_int : [float, float] 46 confidence interval for intercept at confidence level alpha 47 ''' 48 49 # Make sure arrays have the same length 50 assert ( len(X0) == len(Y0) ), 'Arrays X and Y must have the same length' 51 52 # Make sure alpha is within the range 0-1 53 assert (alpha < 1), 'alpha must be less than 1' 54 assert (alpha > 0), 'alpha must be greater than 0' 55 56 # Drop any NaN elements of X or Y 57 # Infinite values are allowed but will make the result undefined 58 idx = ~np.isnan(X0 * Y0 ) 59 X = X0[idx] 60 Y = Y0[idx] 61 62 # Number of finite elements 63 N = len(X) 64 65 # Degrees of freedom for model 66 if intercept: 67 dfmod = 2 68 else: 69 dfmod = 1 70 71 # Robust factors, to be changed later 72 rfac1 = 1 73 rfac2 = 1 74 75 if intercept: 76 # Average values 77 Xmean = np.mean(X) 78 Ymean = np.mean(Y) 79 else: 80 Xmean = 0 81 Ymean = 0 82 83 # Sample variance of X, Y 84 Vx = np.sum( (X - Xmean)**2 ) / (N-1) 85 Vy = np.sum( (Y - Ymean)**2 ) / (N-1) 86 87 # Sample covariance of X, Y 88 Vxy = np.sum( (X - Xmean) * (Y - Ymean) ) / (N-1) 89 90 if robust: 91 92 from sklearn.covariance import MinCovDet 93 94 rcov = MinCovDet().fit( np.array([X,Y]).T ) 95 96 Xmean = rcov.location_[0] 97 Ymean = rcov.location_[1] 98 Vx = rcov.covariance_[0,0] 99 Vy = rcov.covariance_[1,1] 100 Vxy = rcov.covariance_[0,1] 101 102 # Number of observations used in covariance estimate 103 N = rcov.support_.sum() 104 105 if wartonrobust: 106 107 c=np.sqrt(3) 108 109 # Use alternate methods to calculate robust means and covariance 110 q = stats.chi2.cdf(c,2) 111 112 # Huber M to get robust mean and covariance 113 rm, rcov = _huber_cov( np.array([X,Y]), c=c ) 114 115 Xmean = rm[0] 116 Ymean = rm[1] 117 Vx = rcov[0,0] 118 Vy = rcov[1,1] 119 Vxy = rcov[0,1] 120 121 # Robust factors scale the degrees of freedom 122 rfac1, rfac2 = _robust_factor( np.array([X,Y]), rm, rcov, q ) 123 124 # Correlation Coefficient 125 R = Vxy / np.sqrt( Vx * Vy ) 126 127 Sx = np.sqrt( Vx ) 128 Sy = np.sqrt( Vy ) 129 130 # Slope 131 slope = np.sign(R) * Sy / Sx 132 133 # Standard error of slope estimate 134 ste_slope = np.sqrt( rfac2 / (N-dfmod) * Sy**2 / Sx**2 * (1-R**2) ) 135 136 # Confidence interval for slope 137 B = (1-R**2) / (N-dfmod) * rfac1 * stats.f.isf(alpha,1,(N-dfmod)) 138 ci_slope = slope * ( np.sqrt( B+1 ) + np.sqrt(B) * np.array([-1,+1]) ) 139 140 # If slope is negative, flip the order for first element is most negative 141 if slope < 0: 142 ci_slope = np.flipud( ci_slope ) 143 144 if intercept: 145 # Intercept 146 Intercept = Ymean - slope * Xmean 147 148 # Residuals 149 resid = Y - (Intercept + slope * X ) 150 151 # Sample standard deviation of the residuals 152 Sr = np.std( resid, ddof=dfmod ) 153 154 # Another method, may be faster, but less obvious 155 # This method is better for robust fitting, because the simple 156 Sr = np.sqrt((Vy - 2 * slope * Vxy + slope**2 * Vx ) * (N-1) / (N-dfmod) ) 157 #print(Sr,Sr2) 158 159 # Standard error of the intercept estimate 160 ste_int = np.sqrt( Sr**2/N * rfac2 + Xmean**2 * ste_slope**2 ) 161 162 # Confidence interval for Intercept 163 tcrit = stats.t.isf(alpha/2,N-dfmod) 164 ci_int = Intercept + ste_int * np.array([-tcrit,tcrit]) 165 166 else: 167 # Intercept is zero by definition 168 Intercept = 0 169 ste_int = 0 170 ci_int = np.array([0,0]) 171 172 return slope, Intercept, ste_slope, ste_int, ci_slope, ci_int 173 174#@jit 175def _robust_factor( X, rm, rcov, c=0.777 ): 176 '''Robust factor defined by Taskinen and Warton (2013)''' 177 178 rfac1=1 179 rfac2=1 180 181 N = X.shape[1] # number of pointsn 182 k = X.shape[0] # number of variables 183 184 # inverse square root of covariance matrix 185 U, S, V = la.svd( rcov ) 186 rinvsq = U @ np.diag(np.sqrt(1/S)) @ V 187 188 # Means in matrix form 189 Xm = np.tile(rm,[N,1]).T 190 191 # Centered data 192 Xc = X-Xm 193 194 # z-score 195 z = la.norm( rinvsq @ Xc, axis=0 ) 196 197 q = stats.chi2.cdf(3,k) 198 ###CHECK THE C VALUE 199 rfac1 = np.mean( _alpha_fun(z, k, q )**2 ) / 8 200 rfac2 = np.mean( _gamma_fun(z, k, q )**2 ) / 2 201 202 return rfac1, rfac2 203 204def _alpha_fun( r, k, q ): 205 '''Alpha function from Taskinen and Warton (2013)''' 206 207 c = stats.chi2.ppf(q,k) 208 sig = stats.chi2.cdf(c,k+2) + (c/k) * (1-q) 209 c = np.sqrt( c ) 210 211 #c2 = c**2 212 #q = stats.chi2.cdf(c2,k) 213 #s2 = stats.chi2.cdf(c2,k+2) + (c2/k) * (1-q) 214 215 eta = r**2 / (2 * sig**2) 216 eta[r>c] = c**2 / (4 * sig**2) 217 eta = np.mean( eta ) 218 219 alpha = r**2 / (eta * sig**2) 220 alpha[r>c] = c**2 / (eta * sig**2) 221 222 return alpha 223 224def _gamma_fun( r, k, q ): 225 '''Gamma function from Taskinen and Warton (2013)''' 226 227 c = np.sqrt( stats.chi2.ppf(q,k) ) 228 #sig = stats.chi2.cdf(c,k+2) + (c/k) * (1-q) 229 c = np.sqrt( c ) 230 231 eta = c * (k-1) / ( r*k ) 232 eta[r<=c] = 1 233 eta = np.mean(eta) 234 235 gamma = r / eta 236 gamma[r>c] = c/eta 237 238 return gamma 239 240#@jit 241def _huber_cov( X, c=1.73): 242 '''Huber's M robust estimator for mean and covariance 243 244 Method is from Taskinen and Warton, 2013''' 245 246 N = X.shape[1] # number of points 247 k = X.shape[0] # number of variables 248 249 # first guess is normal covariance 250 rcov = np.cov( X ) 251 252 # first guess is normal means 253 rm = X.mean(axis=1) 254 255 # Alternate values 256 #c=3 257 #c=1.345 258 c2 = c**2 259 q = stats.chi2.cdf(c2,k) 260 s2 = stats.chi2.cdf(c2,k+2) + (c2/k) * (1-q) 261 262 # These parameters are used by huber.M, but seem suspicious 263 c=3 264 q = stats.chi2.cdf(c,k) 265 c2=stats.chi2.cdf(c,k+2) + (c/k) * (1-q) 266 267 R = la.cholesky( la.inv( rcov ) ) 268 269 #q = 0.777, c = 3 270 271 it = 0 272 d1 = 1 273 d2 = 1 274 eps = 1e-6 275 maxit = 100 276 while ((it<maxit) and (d1>eps) and (d2>eps) ): 277 278 # Means in matrix form 279 Xm = np.tile(rm,[N,1]).T 280 281 # Centered X 282 Xc = X - Xm 283 284 # z scores 285 s = np.diag( Xc.T @ (R.T @ R ) @ Xc ) 286 287 # Weighting 288 u = (c/c2) / s 289 u[s<=c] = 1/c2 290 291 # updated weighted inverse covariance 292 C = R @ (Xc @ np.diag(u) @ Xc.T) @ R.T / N 293 294 R0 = la.cholesky( la.inv( C ) ) 295 296 R = R0 @ R 297 298 d1 = np.max( np.abs( R0 - np.identity(k)).sum(axis=1) ) 299 300 # updated z scores 301 s = np.diag( Xc.T @ (R.T @ R) @ Xc ) 302 303 # Weighting 304 v = np.sqrt(c / s ) 305 v[s<=c] = 1 306 307 # Update increment for means 308 h = ( Xc @ np.diag(v) ).mean(axis=1) / v.mean() 309 310 # Updated robust mean 311 rm = rm + h 312 313 d2 = np.sqrt( (h-rm).T @ (h-rm) ) 314 315 it += 1 316 317 if it>=maxit: 318 raise RuntimeError('huber_cov did not converge') 319 320 # Robust covariance 321 rcov = la.inv( R.T @ R ) 322 323 return rm, rcov 324 325def _huber_cov_old( X, c=1.73): 326 '''Huber's M robust estimator for mean and covariance, old version 327 328 Method is from Taskinen and Warton, 2013''' 329 330 N = X.shape[1] # number of points 331 k = X.shape[0] # number of variables 332 333 # first guess is normal covariance 334 rcov = np.cov( X ) 335 336 # first guess is normal means 337 rm = X.mean(axis=1) 338 339 # Alternate values 340 #c=3 341 #c=1.345 342 c2 = c**2 343 q = stats.chi2.cdf(c2,k) 344 s2 = stats.chi2.cdf(c2,k+2) + (c2/k) * (1-q) 345 346 #q = 0.777, c = 3 347 348 it = 0 349 d1 = 1 350 d2 = 1 351 eps = 1e-6 352 while ((it<100) and (d1>eps) and (d2>eps) ): 353 354 # inverse square root of covariance matrix 355 U, S, V = la.svd( rcov ) 356 rinvsq = U @ np.diag(np.sqrt(1/S)) @ V 357 #rinvsq = la.sqrtm( la.inv( rcov ) ) 358 359 # z-scores 360 z = la.norm( rinvsq @ ( X - np.tile(rm,[N,1]).T ) , axis=0 ) 361 362 # Weights for mean 363 w1 = np.minimum( 1, c/z ) 364 #w1 = np.ones(N) 365 366 # Weights for covariance, add scaling later 367 w2 = w1**2 / s2 368 369 # New estimate of robust mean 370 rm2 = np.sum( X * np.tile(w1,[2,1]), axis=1 ) / np.sum(w1) 371 372 # Means in matrix form 373 Xm = np.tile(rm2,[N,1]).T 374 375 # New estimate of robust covariance 376 rcov2 = 1/(N-1) * ( (X-Xm) @ ( (X-Xm) * np.tile(w2,[2,1]) ).T ) 377 378 it += 1 379 d1 = np.max(np.abs(rm2-rm)) 380 d2 = np.max(np.abs(rcov2-rcov)) 381 rcov = rcov2 382 rm = rm2 383 384 if it>100: 385 raise RuntimeError('huber_cov_old did not converge') 386 387 return rm, rcov
def
sma_warton(X0, Y0, alpha=0.05, intercept=True, robust=False, wartonrobust=False):
20def sma_warton(X0,Y0,alpha=0.05,intercept=True,robust=False,wartonrobust=False): 21 '''Standard major axis, aka reduced major axis 22 23 The main advantage of this is that the best fit 24 of Y to X will be the same as the best fit of X to Y. 25 More details in Warton et al. Biology Review 2006 26 27 Parameters 28 ---------- 29 X, Y : array_like 30 Input values, Must have same length. 31 alpha : float 32 Desired confidence level for output. 33 34 Returns 35 ------- 36 Gradient : float 37 Gradient or Slope of Y vs. X 38 Intercept : floats 39 Y intercept. 40 stdgrad : float 41 Standard error of gradient estimate 42 stdint : float 43 standard error of intercept estimate 44 ci_grad : [float, float] 45 confidence interval for gradient at confidence level alpha 46 ci_int : [float, float] 47 confidence interval for intercept at confidence level alpha 48 ''' 49 50 # Make sure arrays have the same length 51 assert ( len(X0) == len(Y0) ), 'Arrays X and Y must have the same length' 52 53 # Make sure alpha is within the range 0-1 54 assert (alpha < 1), 'alpha must be less than 1' 55 assert (alpha > 0), 'alpha must be greater than 0' 56 57 # Drop any NaN elements of X or Y 58 # Infinite values are allowed but will make the result undefined 59 idx = ~np.isnan(X0 * Y0 ) 60 X = X0[idx] 61 Y = Y0[idx] 62 63 # Number of finite elements 64 N = len(X) 65 66 # Degrees of freedom for model 67 if intercept: 68 dfmod = 2 69 else: 70 dfmod = 1 71 72 # Robust factors, to be changed later 73 rfac1 = 1 74 rfac2 = 1 75 76 if intercept: 77 # Average values 78 Xmean = np.mean(X) 79 Ymean = np.mean(Y) 80 else: 81 Xmean = 0 82 Ymean = 0 83 84 # Sample variance of X, Y 85 Vx = np.sum( (X - Xmean)**2 ) / (N-1) 86 Vy = np.sum( (Y - Ymean)**2 ) / (N-1) 87 88 # Sample covariance of X, Y 89 Vxy = np.sum( (X - Xmean) * (Y - Ymean) ) / (N-1) 90 91 if robust: 92 93 from sklearn.covariance import MinCovDet 94 95 rcov = MinCovDet().fit( np.array([X,Y]).T ) 96 97 Xmean = rcov.location_[0] 98 Ymean = rcov.location_[1] 99 Vx = rcov.covariance_[0,0] 100 Vy = rcov.covariance_[1,1] 101 Vxy = rcov.covariance_[0,1] 102 103 # Number of observations used in covariance estimate 104 N = rcov.support_.sum() 105 106 if wartonrobust: 107 108 c=np.sqrt(3) 109 110 # Use alternate methods to calculate robust means and covariance 111 q = stats.chi2.cdf(c,2) 112 113 # Huber M to get robust mean and covariance 114 rm, rcov = _huber_cov( np.array([X,Y]), c=c ) 115 116 Xmean = rm[0] 117 Ymean = rm[1] 118 Vx = rcov[0,0] 119 Vy = rcov[1,1] 120 Vxy = rcov[0,1] 121 122 # Robust factors scale the degrees of freedom 123 rfac1, rfac2 = _robust_factor( np.array([X,Y]), rm, rcov, q ) 124 125 # Correlation Coefficient 126 R = Vxy / np.sqrt( Vx * Vy ) 127 128 Sx = np.sqrt( Vx ) 129 Sy = np.sqrt( Vy ) 130 131 # Slope 132 slope = np.sign(R) * Sy / Sx 133 134 # Standard error of slope estimate 135 ste_slope = np.sqrt( rfac2 / (N-dfmod) * Sy**2 / Sx**2 * (1-R**2) ) 136 137 # Confidence interval for slope 138 B = (1-R**2) / (N-dfmod) * rfac1 * stats.f.isf(alpha,1,(N-dfmod)) 139 ci_slope = slope * ( np.sqrt( B+1 ) + np.sqrt(B) * np.array([-1,+1]) ) 140 141 # If slope is negative, flip the order for first element is most negative 142 if slope < 0: 143 ci_slope = np.flipud( ci_slope ) 144 145 if intercept: 146 # Intercept 147 Intercept = Ymean - slope * Xmean 148 149 # Residuals 150 resid = Y - (Intercept + slope * X ) 151 152 # Sample standard deviation of the residuals 153 Sr = np.std( resid, ddof=dfmod ) 154 155 # Another method, may be faster, but less obvious 156 # This method is better for robust fitting, because the simple 157 Sr = np.sqrt((Vy - 2 * slope * Vxy + slope**2 * Vx ) * (N-1) / (N-dfmod) ) 158 #print(Sr,Sr2) 159 160 # Standard error of the intercept estimate 161 ste_int = np.sqrt( Sr**2/N * rfac2 + Xmean**2 * ste_slope**2 ) 162 163 # Confidence interval for Intercept 164 tcrit = stats.t.isf(alpha/2,N-dfmod) 165 ci_int = Intercept + ste_int * np.array([-tcrit,tcrit]) 166 167 else: 168 # Intercept is zero by definition 169 Intercept = 0 170 ste_int = 0 171 ci_int = np.array([0,0]) 172 173 return slope, Intercept, ste_slope, ste_int, ci_slope, ci_int
Standard major axis, aka reduced major axis
The main advantage of this is that the best fit of Y to X will be the same as the best fit of X to Y. More details in Warton et al. Biology Review 2006
Parameters
- X, Y (array_like): Input values, Must have same length.
- alpha (float): Desired confidence level for output.
Returns
- Gradient (float): Gradient or Slope of Y vs. X
- Intercept (floats): Y intercept.
- stdgrad (float): Standard error of gradient estimate
- stdint (float): standard error of intercept estimate
- ci_grad ([float, float]): confidence interval for gradient at confidence level alpha
- ci_int ([float, float]): confidence interval for intercept at confidence level alpha