scipy.stats.

chatterjeexi#

scipy.stats.chatterjeexi(x, y, *, axis=0, y_continuous=False, method='asymptotic', nan_policy='propagate', keepdims=False)[source]#

Compute the xi correlation and perform a test of independence

The xi correlation coefficient is a measure of association between two variables; the value tends to be close to zero when the variables are independent and close to 1 when there is a strong association. Unlike other correlation coefficients, the xi correlation is effective even when the association is not monotonic.

Parameters:
x, yarray-like

The samples: corresponding observations of the independent and dependent variable. The (N-d) arrays must be broadcastable.

axisint or None, default: 0

If an int, the axis of the input along which to compute the statistic. The statistic of each axis-slice (e.g. row) of the input will appear in a corresponding element of the output. If None, the input will be raveled before computing the statistic.

method‘asymptotic’ or PermutationMethod instance, optional

Selects the method used to calculate the p-value. Default is ‘asymptotic’. The following options are available.

  • 'asymptotic': compares the standardized test statistic against the normal distribution.

  • PermutationMethod instance. In this case, the p-value is computed using permutation_test with the provided configuration options and other appropriate settings.

y_continuousbool, default: False

Whether y is assumed to be drawn from a continuous distribution. If y is drawn from a continuous distribution, results are valid whether this is assumed or not, but enabling this assumption will result in faster computation and typically produce similar results.

nan_policy{‘propagate’, ‘omit’, ‘raise’}

Defines how to handle input NaNs.

  • propagate: if a NaN is present in the axis slice (e.g. row) along which the statistic is computed, the corresponding entry of the output will be NaN.

  • omit: NaNs will be omitted when performing the calculation. If insufficient data remains in the axis slice along which the statistic is computed, the corresponding entry of the output will be NaN.

  • raise: if a NaN is present, a ValueError will be raised.

keepdimsbool, default: False

If this is set to True, the axes which are reduced are left in the result as dimensions with size one. With this option, the result will broadcast correctly against the input array.

Returns:
resSignificanceResult

An object containing attributes:

statisticfloat

The xi correlation statistic.

pvaluefloat

The associated p-value: the probability of a statistic at least as high as the observed value under the null hypothesis of independence.

Notes

There is currently no special handling of ties in x; they are broken arbitrarily by the implementation. [1] recommends: “if there are ties among the Xi’s, then choose an increasing rearrangement as above by breaking ties uniformly at random.” This is easily accomplished by adding a small amount of random noise to x; see examples.

[1] notes that the statistic is not symmetric in x and y by design: “…we may want to understand if \(Y\) is a function \(X\), and not just if one of the variables is a function of the other.” See [1] Remark 1.

Beginning in SciPy 1.9, np.matrix inputs (not recommended for new code) are converted to np.ndarray before the calculation is performed. In this case, the output will be a scalar or np.ndarray of appropriate shape rather than a 2D np.matrix. Similarly, while masked elements of masked arrays are ignored, the output will be a scalar or np.ndarray rather than a masked array with mask=False.

Array API Standard Support

chatterjeexi has experimental support for Python Array API Standard compatible backends in addition to NumPy. Please consider testing these features by setting an environment variable SCIPY_ARRAY_API=1 and providing CuPy, PyTorch, JAX, or Dask arrays as array arguments. The following combinations of backend and device (or other capability) are supported.

Library

CPU

GPU

NumPy

n/a

CuPy

n/a

PyTorch

JAX

Dask

n/a

See Support for the array API standard for more information.

References

[1] (1,2,3,4,5)

Chatterjee, Sourav. “A new coefficient of correlation.” Journal of the American Statistical Association 116.536 (2021): 2009-2022. DOI:10.1080/01621459.2020.1758115.

Examples

Generate perfectly correlated data, and observe that the xi correlation is nearly 1.0.

>>> import numpy as np
>>> from scipy import stats
>>> rng = np.random.default_rng()
>>> x = rng.uniform(0, 10, size=100)
>>> y = np.sin(x)
>>> res = stats.chatterjeexi(x, y)
>>> res.statistic
np.float64(0.9012901290129013)

The probability of observing such a high value of the statistic under the null hypothesis of independence is very low.

>>> res.pvalue
np.float64(2.2206974648177804e-46)

As noise is introduced, the correlation coefficient decreases.

>>> noise = rng.normal(scale=[[0.1], [0.5], [1]], size=(3, 100))
>>> res = stats.chatterjeexi(x, y + noise, axis=-1)
>>> res.statistic
array([0.79507951, 0.41824182, 0.16651665])

Because the distribution of y is continuous, it is valid to pass y_continuous=True. The statistic is identical, and the p-value (not shown) is only slightly different.

>>> stats.chatterjeexi(x, y + noise, y_continuous=True, axis=-1).statistic
array([0.79507951, 0.41824182, 0.16651665])

Consider a case in which there are ties in x.

>>> x = rng.integers(10, size=1000)
>>> y = rng.integers(10, size=1000)

[1] recommends breaking the ties uniformly at random.

>>> d = rng.uniform(1e-5, size=x.size)
>>> res = stats.chatterjeexi(x + d, y)
>>> res.statistic
-0.029919991638798438

Since this gives a randomized estimate of the statistic, [1] also suggests considering the average over all possibilities of breaking ties. This is computationally infeasible when there are many ties, but a randomized estimate of this quantity can be obtained by considering many random possibilities of breaking ties.

>>> d = rng.uniform(1e-5, size=(9999, x.size))
>>> res = stats.chatterjeexi(x + d, y, axis=1)
>>> np.mean(res.statistic)
0.001186895213756626