corr_sig (pyleoclim.utils.correlation.corr_sig)¶
- pyleoclim.utils.correlation.corr_sig(y1, y2, nsim=1000, method='isospectral', alpha=0.05)[source]¶
Estimates the Pearson’s correlation and associated significance between two non IID time series
The significance of the correlation is assessed using one of the following methods:
- ‘ttest’: T-test adjusted for effective sample size.
This is a parametric test (data are Gaussian and identically distributed) with a rather ad-hoc adjustment. It is instantaneous but makes a lot of assumptions about the data, many of which may not be met.
- ‘isopersistent’: AR(1) modeling of x and y.
This is a parametric test as well (series follow an AR(1) model) but solves the issue by direct simulation.
- ‘isospectral’: phase randomization of original inputs. (default)
This is a non-parametric method, assuming only wide-sense stationarity.
For 2 and 3, computational requirements scale with nsim. When possible, nsim should be at least 1000.
- Parameters
y1 (array) – vector of (real) numbers of same length as y2, no NaNs allowed
y2 (array) – vector of (real) numbers of same length as y1, no NaNs allowed
nsim (int) – the number of simulations [default: 1000]
method ({'ttest','isopersistent','isospectral' (default)}) – method for significance testing
alpha (float) – significance level for critical value estimation [default: 0.05]
- Returns
res – the result dictionary, containing
- rfloat
correlation coefficient
- pfloat
the p-value
- signifbool
true if significant; false otherwise Note that signif = True if and only if p <= alpha.
- Return type
dict
See also
pyleoclim.utils.correlation.corr_ttest
Estimates the significance of correlations between 2 time series using the classical T-test adjusted for effective sample size.
pyleoclim.utils.correlation.corr_isopersist
Computes correlation between two timeseries, and their significance using Ar(1) modeling.
pyleoclim.utils.correlation.corr_isospec
Estimates the significance of the correlation using phase randomization