correlation_fft

mdhelper.algorithm.correlation.correlation_fft(arr1: ndarray[float], arr2: ndarray[float] = None, axis: int = None, *, average: bool = False, double: bool = False, vector: bool = False) ndarray[float][source]

Evaluates the autocorrelation function (ACF) or cross-correlation function (CCF) of a time series using fast Fourier transforms (FFT).

The algorithm, better known as the Fast Correlation Algorithm (FCA) [1] [2], is a result of the Wiener–Khinchin theorem and has a time complexity of O(NlogN). Effectively, the ACF can be computed from the raw data r(t) with two FFTs using

r^(ξ)=FFT(r(t))A(τ)=FFT1(r^(ξ)r^(ξ))

where τ is the time lag and the asterisk () denotes the complex conjugate.

Similarly, the CCF for species i and j is evaluated using

Cij(τ)=FFT1(FFT(ri(t))FFT(rj(t)))
Parameters:
arr1numpy.ndarray

Time evolution of N entities over Nb blocks of Nt frames each.

Shape:

  • Scalar: (Nt,), (Nt,N), (Nb,Nt), or (Nb,Nt,N).

  • Vector: (Nt,d), (Nt,N,Nd), (Nb,Nt,Nd), or (Nb,Nt,N,Nd), where Nd is the number of dimensions each vector has.

arr2numpy.ndarray, optional

Time evolution of another N entities. If provided, the CCF for arr1 and arr2 is calculated. Otherwise, the ACF for arr1 is calculated.

Shape: Same as arr1.

axisint, optional

Axis along which to evaluate the ACF/CCF. If arr1 contains a full, unsplit trajectory, the ACF/CCF should be evaluated along the first axis (axis=0). If arr1 contains a trajectory split into multiple blocks, the ACF/CCF should be evaluated along the second axis (axis=1). If not specified, the axis is determined automatically using the shape of arr1.

averagebool, keyword-only, default: True

Determines whether the ACF/CCF is averaged over all entities if the arrays contain information for multiple entities.

doublebool, keyword-only, default: False

If True, the ACF is doubled or the CCFs for the negative and positive time lags are combined. Useful for evaluating the mean squared or cross displacement. See mdhelper.algorithm.correlation.msd_fft() for more information.

vectorbool, keyword-only, default: False

Specifies whether arr1 and arr2 contain vectors. If True, the ACF/CCF is summed over the last dimension.

Returns:
corrnumpy.ndarray

Autocorrelation or cross-correlation function.

Shape:

For ACF, the shape is that of arr1 but with the following modifications:

  • If average=True, the axis containing the N entities is removed.

  • If vector=True, the last dimension is removed.

For CCF, the shape is that of arr1 but with the following modifications:

  • If average=True, the axis containing the N entities is removed.

  • If double=False, the axis containing the Nt times now has a length of 2Nt1 to accomodate negative and positive time lags.

  • If vector=True, the last dimension is removed.

References

[1]

Kneller, G. R.; Keiner, V.; Kneller, M.; Schiller, M. NMOLDYN: A Program Package for a Neutron Scattering Oriented Analysis of Molecular Dynamics Simulations. Computer Physics Communications 1995, 91 (1–3), 191–214. https://doi.org/10.1016/0010-4655(95)00048-K.

[2]

Calandrini, V.; Pellegrini, E.; Calligari, P.; Hinsen, K.; Kneller, G. R. NMoldyn - Interfacing Spectroscopic Experiments, Molecular Dynamics Simulations and Models for Time Correlation Functions. JDN 2011, 12, 201–232. https://doi.org/10.1051/sfn/201112010.