spatialstats.particles

Spatial distribution functions

Routines to calculate the 2-point spatial distribution function \(g(\mathbf{r})\) for displacements between point and rod-like particle pairs. Can reduce to the usual isotropic \(g(r)\) with the corresponding structure factor \(S(q)\).

See here to learn more.

Adapted from https://github.com/wenyan4work/point_cloud.

sdf(positions, boxsize, orientations=None, rmin=None, rmax=None, nr=100, nphi=None, ntheta=None, int=<class 'numpy.int32'>, float=<class 'numpy.float64'>, bench=False, **kwargs)

Calculate the spatial distribution function \(g(r, \phi)\) or \(g(r, \phi, \theta)\) for a set of \(N\) point-like or rod-like particles \(\mathbf{r}_i\) in a 2D or 3D periodic box. \((r, \phi, \theta)\) are the spherical coordinates for displacement vectors between particle pairs \(\mathbf{r}_{ij} = \mathbf{r}_i - \mathbf{r}_j\). \(\phi\) is the azimuthal angle and \(\theta\) is the inclination angle.

If particles orientations \(\mathbf{p}_i\) are included, instead define \((r, \phi, \theta)\) as the coordinate system with \(\mathbf{p}_i\) pointed in the \(+z\) direction.

Note

Reduces to the 1D distribution function \(g(r)\) when nphi = None and ntheta = None.

Parameters
  • positions (np.ndarray, shape (N, ndim)) – Particle positions \(\mathbf{r}_i\) in 2D or 3D for \(N\) particles. Passed to scipy.spatial.cKDTree.

  • boxsize (list of float) – The rectangular domain over which to apply periodic boundary conditions. Passed to scipy.spatial.cKDTree.

  • orientations (np.ndarray, shape (N, ndim), optional) – Particle orientation vectors \(\mathbf{p}_i\). Vectors should be unitary, but they will be normalized automatically.

  • rmin (float, optional) – Minimum \(r\) value in \(g(r, \phi, \theta)\).

  • rmax (float, optional) – Cutoff radius for KDTree search and maximum \(r\) value in \(g(r, \phi, \theta)\). Default is half the maximum dimension of boxsize.

  • nr (int, optional) – Number of points to bin in \(r\).

  • nphi (int, optional) – Number of points to bin in \(\phi\).

  • ntheta (int, optional) – Number of points to bin in \(\theta\).

  • int (np.dtype, optional) – Integer type for pair counting array. Lets the user relax memory requirements.

  • float (np.dtype, optional) – Floating-point type for displacement buffers. Lets the user relax memory requirements.

  • bench (bool, optional) – Print message for time of calculation.

Returns

  • g (np.ndarray, shape (nr, nphi, ntheta)) – Radial distribution function \(g(r, \phi, \theta)\). If the user does not bin for a certain coordinate, g will not be 3 dimensional (e.g. if nphi = None, g will be shape (nr, ntheta)).

  • r (np.ndarray, shape (nr,)) – Left edges of radial bins \(r\).

  • phi (np.ndarray, shape (nphi,)) – Left edges of angular bins \(\phi \in [-\pi, \pi)\).

  • theta (np.ndarray, shape (ntheta,)) – Left edges of angular bins \(\theta \in [0, \pi)\). Not returned for 2D datasets.

structure_factor(gr, r, N, boxsize, q=None, **kwargs)

Calculate the isotropic structure factor \(S(q)\) from the radial distribution function \(g(r)\) of a set of \(N\) particle positions in a 2D or 3D periodic box with volume \(V\).

The structure factor in 3D is computed as

\[S(q) = 1 + 4\pi \rho \frac{1}{q} \int dr \ r \ \textrm{sin}(qr) [g(r) - 1]\]

and in 2D

\[S(q) = 1 + 2\pi \rho \int dr \ r \ J_{0}(qr) [g(r) - 1],\]

where \(\rho = N/V\) and \(J_{0}\) is the 0th bessel function of the first kind.

Parameters
  • gr (np.ndarray) – The radial distribution function \(g(r)\) from spatialstats.particles.sdf.

  • r (np.ndarray) – The domain of \(g(r)\) from spatialstats.particles.sdf.

  • N (int) – The number of particles \(N\).

  • boxsize (list of float) – The rectangular domain over which to apply periodic boundary conditions. See scipy.spatial.cKDTree.

  • q (np.ndarray, optional) – Dimensional wavenumber bins \(q\). If None, q = np.arange(dq, 200*dq, dq) with dq = 2*np.pi / max(boxsize).

Returns

  • Sq (np.ndarray) – The structure factor \(S(q)\).

  • q (np.ndarray) – Wavenumber bins \(q\).