MSD Module¶
Overview
Compute the mean squared displacement. |
Details
The freud.msd
module provides functions for computing the
mean-squared-displacement (MSD) of particles in periodic systems.
- class freud.msd.MSD¶
Bases:
_Compute
Compute the mean squared displacement.
The mean squared displacement (MSD) measures how much particles move over time. The MSD plays an important role in characterizing Brownian motion, since it provides a measure of whether particles are moving according to diffusion alone or if there are other forces contributing. There are a number of definitions for the mean squared displacement. This function provides access to the two most common definitions through the mode argument.
'window'
(default): This mode calculates the most common form of the MSD, which is defined as\[MSD(m) = \frac{1}{N_{particles}} \sum_{i=1}^{N_{particles}} \frac{1}{N-m} \sum_{k=0}^{N-m-1} (\vec{r}_i(k+m) - \vec{r}_i(k))^2\]where \(r_i(t)\) is the position of particle \(i\) in frame \(t\). According to this definition, the mean squared displacement is the average displacement over all windows of length \(m\) over the course of the simulation. Therefore, for any \(m\), \(MSD(m)\) is averaged over all windows of length \(m\) and over all particles. This calculation can be accessed using the ‘window’ mode of this function.
The windowed calculation can be quite computationally intensive. To perform this calculation efficiently, we use the algorithm described in [CPC+11] as described in this StackOverflow thread.
Note
The most intensive part of this calculation is computing an FFT. To maximize performance, freud attempts to use the fastest FFT library available. By default, the order of preference is pyFFTW, SciPy, and then NumPy. If you are experiencing significant slowdowns in calculating the MSD, you may benefit from installing a faster FFT library, which freud will automatically detect. The performance change will be especially noticeable if the length of your trajectory is a number whose prime factorization consists of extremely large prime factors. The standard Cooley-Tukey FFT algorithm performs very poorly in this case, so installing pyFFTW will significantly improve performance.
Note that while pyFFTW is released under the BSD 3-Clause license, the FFTW library is available under either GPL or a commercial license.
'direct'
: Under some circumstances, however, we may be more interested in calculating a different quantity described by\begin{eqnarray*} MSD(t) =& \dfrac{1}{N_{particles}} \sum_{i=1}^{N_{particles}} (r_i(t) - r_i(0))^2 \\ \end{eqnarray*}In this case, at each time point (i.e. simulation frame) we simply compute how much particles have moved from their initial position, averaged over all particles. For more information on this calculation, see the Wikipedia page.
Note
The MSD is only well-defined when the box is constant over the course of the simulation. Additionally, the number of particles must be constant over the course of the simulation.
- Parameters:
box (
freud.box.Box
, optional) – If not provided, the class will assume that all positions provided in calls tocompute()
are already unwrapped. (Default value =None
).mode (str, optional) – Mode of calculation. Options are
'window'
and'direct'
. (Default value ='window'
).
- box¶
Box used in the calculation.
- Type:
- compute(self, positions, images=None, reset=True)¶
Calculate the MSD for the positions provided.
Note
Unlike most methods in freud, accumulation for the MSD is split over points rather than frames of a simulation. The reason for this choice is that efficient computation of the MSD requires using the entire trajectory for a given particle. As a result, when setting
reset=False
, you must provide the positions of each point over the full length of the trajectory, but you may callcompute
multiple times with different subsets the points to calculate the MSD over the full set of positions. The primary use-case is when the trajectory is so large that computing an MSD on all particles at once is prohibitively expensive.- Parameters:
positions ((\(N_{frames}\), \(N_{particles}\), 3)
numpy.ndarray
) – The particle positions over a trajectory. If neither box nor images are provided, the positions are assumed to be unwrapped already.images ((\(N_{frames}\), \(N_{particles}\), 3)
numpy.ndarray
, optional) – The particle images to unwrap with if provided. Must be provided along with a simulation box (in the constructor) if particle positions need to be unwrapped. If neither are provided, positions are assumed to be unwrapped already. (Default value =None
).reset (bool) – Whether to erase the previously computed values before adding the new computation; if False, will accumulate data (Default value: True).
- property msd¶
The mean squared displacement.
- Type:
\(\left(N_{frames}, \right)\)
numpy.ndarray
- property particle_msd¶
The per particle based mean squared displacement.
- Type:
\(\left(N_{frames}, N_{particles} \right)\)
numpy.ndarray
- plot(self, ax=None)¶
Plot MSD.
- Parameters:
ax (
matplotlib.axes.Axes
, optional) – Axis to plot on. IfNone
, make a new figure and axis. (Default value =None
)- Returns:
Axis with the plot.
- Return type: