Order Module

Overview

freud.order.CubaticOrderParameter

Compute the cubatic order parameter [HajiAkbari2015] for a system of particles using simulated annealing instead of Newton-Raphson root finding.

freud.order.NematicOrderParameter

Compute the nematic order parameter for a system of particles.

freud.order.HexOrderParameter

Calculates the \(k\)-atic order parameter for each particle in the system.

freud.order.TransOrderParameter

Compute the translational order parameter for each particle.

freud.order.LocalQl

Compute the local Steinhardt [Steinhardt1983] rotationally invariant \(Q_l\) order parameter for a set of points.

freud.order.LocalQlNear

A variant of the LocalQl class that performs its average over nearest neighbor particles as determined by an instance of freud.locality.NeighborList.

freud.order.LocalWl

Compute the local Steinhardt [Steinhardt1983] rotationally invariant \(W_l\) order parameter for a set of points.

freud.order.LocalWlNear

A variant of the LocalWl class that performs its average over nearest neighbor particles as determined by an instance of freud.locality.NeighborList.

freud.order.SolLiq

Uses dot products of \(Q_{lm}\) between particles for clustering.

freud.order.SolLiqNear

A variant of the SolLiq class that performs its average over nearest neighbor particles as determined by an instance of freud.locality.NeighborList.

freud.order.RotationalAutocorrelation

Calculates a measure of total rotational autocorrelation based on hyperspherical harmonics as laid out in “Design rules for engineering colloidal plastic crystals of hard polyhedra - phase behavior and directional entropic forces” by Karas et al.

Details

The freud.order module contains functions which compute order parameters for the whole system or individual particles. Order parameters take bond order data and interpret it in some way to quantify the degree of order in a system using a scalar value. This is often done through computing spherical harmonics of the bond order diagram, which are the spherical analogue of Fourier Transforms.

Cubatic Order Parameter

class freud.order.CubaticOrderParameter(t_initial, t_final, scale, n_replicates, seed)

Compute the cubatic order parameter [HajiAkbari2015] for a system of particles using simulated annealing instead of Newton-Raphson root finding.

Module author: Eric Harper <harperic@umich.edu>

Parameters
  • t_initial (float) – Starting temperature.

  • t_final (float) – Final temperature.

  • scale (float) – Scaling factor to reduce temperature.

  • n_replicates (unsigned int) – Number of replicate simulated annealing runs.

  • seed (unsigned int) – Random seed to use in calculations. If None, system time is used.

Variables
  • t_initial (float) – The value of the initial temperature.

  • t_final (float) – The value of the final temperature.

  • scale (float) – The scale

  • cubatic_order_parameter (float) – The cubatic order parameter.

  • orientation (\(\left(4 \right)\) numpy.ndarray) – The quaternion of global orientation.

  • particle_order_parameter (numpy.ndarray) – Cubatic order parameter.

  • particle_tensor (\(\left(N_{particles}, 3, 3, 3, 3 \right)\) numpy.ndarray) – Rank 5 tensor corresponding to each individual particle orientation.

  • global_tensor (\(\left(3, 3, 3, 3 \right)\) numpy.ndarray) – Rank 4 tensor corresponding to global orientation.

  • cubatic_tensor (\(\left(3, 3, 3, 3 \right)\) numpy.ndarray) – Rank 4 cubatic tensor.

  • gen_r4_tensor (\(\left(3, 3, 3, 3 \right)\) numpy.ndarray) – Rank 4 tensor corresponding to each individual particle orientation.

compute

Calculates the per-particle and global order parameter.

Parameters

orientations ((\(N_{particles}\), 4) numpy.ndarray) – Orientations as angles to use in computation.

Nematic Order Parameter

class freud.order.NematicOrderParameter(u)

Compute the nematic order parameter for a system of particles.

Module author: Jens Glaser <jsglaser@umich.edu>

New in version 0.7.0.

Parameters

u (\(\left(3 \right)\) numpy.ndarray) – The nematic director of a single particle in the reference state (without any rotation applied).

Variables
  • nematic_order_parameter (float) – Nematic order parameter.

  • director (\(\left(3 \right)\) numpy.ndarray) – The average nematic director.

  • particle_tensor (\(\left(N_{particles}, 3, 3 \right)\) numpy.ndarray) – One 3x3 matrix per-particle corresponding to each individual particle orientation.

  • nematic_tensor (\(\left(3, 3 \right)\) numpy.ndarray) – 3x3 matrix corresponding to the average particle orientation.

compute

Calculates the per-particle and global order parameter.

Parameters

orientations (\(\left(N_{particles}, 4 \right)\) numpy.ndarray) – Orientations to calculate the order parameter.

Hexatic Order Parameter

class freud.order.HexOrderParameter(rmax, k, n)

Calculates the \(k\)-atic order parameter for each particle in the system.

The \(k\)-atic order parameter for a particle \(i\) and its \(n\) neighbors \(j\) is given by:

\(\psi_k \left( i \right) = \frac{1}{n} \sum_j^n e^{k i \phi_{ij}}\)

The parameter \(k\) governs the symmetry of the order parameter while the parameter \(n\) governs the number of neighbors of particle \(i\) to average over. \(\phi_{ij}\) is the angle between the vector \(r_{ij}\) and \(\left( 1,0 \right)\).

Note

2D: freud.order.HexOrderParameter properly handles 2D boxes. The points must be passed in as [x, y, 0]. Failing to set z=0 will lead to undefined behavior.

Module author: Eric Harper <harperic@umich.edu>

Parameters
  • rmax (float) – +/- r distance to search for neighbors.

  • k (unsigned int) – Symmetry of order parameter (\(k=6\) is hexatic).

  • n (unsigned int) – Number of neighbors (\(n=k\) if \(n\) not specified).

Variables
  • psi (\(\left(N_{particles} \right)\) numpy.ndarray) – Order parameter.

  • box (freud.box.Box) – Box used in the calculation.

  • num_particles (unsigned int) – Number of particles.

  • K (unsigned int) – Symmetry of the order parameter.

compute

Calculates the correlation function and adds to the current histogram.

Parameters

Translational Order Parameter

class freud.order.TransOrderParameter(rmax, k, n)

Compute the translational order parameter for each particle.

Module author: Wenbo Shen <shenwb@umich.edu>

Parameters
  • rmax (float) – +/- r distance to search for neighbors.

  • k (float) – Symmetry of order parameter (\(k=6\) is hexatic).

  • n (unsigned int) – Number of neighbors (\(n=k\) if \(n\) not specified).

Variables
  • d_r (\(\left(N_{particles}\right)\) numpy.ndarray) – Reference to the last computed translational order array.

  • box (freud.box.Box) – Box used in the calculation.

  • num_particles (unsigned int) – Number of particles.

  • K (float) – Normalization value (d_r is divided by K).

compute

Calculates the local descriptors.

Parameters

Steinhardt \(Q_l\) Order Parameter

class freud.order.LocalQl(box, rmax, l, rmin)

Compute the local Steinhardt [Steinhardt1983] rotationally invariant \(Q_l\) order parameter for a set of points.

Implements the local rotationally invariant \(Q_l\) order parameter described by Steinhardt. For a particle i, we calculate the average \(Q_l\) by summing the spherical harmonics between particle \(i\) and its neighbors \(j\) in a local region: \(\overline{Q}_{lm}(i) = \frac{1}{N_b} \displaystyle\sum_{j=1}^{N_b} Y_{lm}(\theta(\vec{r}_{ij}), \phi(\vec{r}_{ij}))\). The particles included in the sum are determined by the rmax argument to the constructor.

This is then combined in a rotationally invariant fashion to remove local orientational order as follows: \(Q_l(i)=\sqrt{\frac{4\pi}{2l+1} \displaystyle\sum_{m=-l}^{l} |\overline{Q}_{lm}|^2 }\).

The computeAve() method provides access to a variant of this parameter that performs a average over the first and second shell combined [Lechner2008]. To compute this parameter, we perform a second averaging over the first neighbor shell of the particle to implicitly include information about the second neighbor shell. This averaging is performed by replacing the value \(\overline{Q}_{lm}(i)\) in the original definition by the average value of \(\overline{Q}_{lm}(k)\) over all the \(k\) neighbors of particle \(i\) as well as itself.

The computeNorm() and computeAveNorm() methods provide normalized versions of compute() and computeAve(), where the normalization is performed by dividing by the average \(Q_{lm}\) values over all particles.

Module author: Xiyu Du <xiyudu@umich.edu>

Module author: Vyas Ramasubramani <vramasub@umich.edu>

Parameters
  • box (freud.box.Box) – Simulation box.

  • rmax (float) – Cutoff radius for the local order parameter. Values near the first minimum of the RDF are recommended.

  • l (unsigned int) – Spherical harmonic quantum number l. Must be a positive integer.

  • rmin (float) – Lower bound for computing the local order parameter. Allows looking at, for instance, only the second shell, or some other arbitrary RDF region (Default value = 0).

Variables
  • box (freud.box.Box) – Box used in the calculation.

  • num_particles (unsigned int) – Number of particles.

  • Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(Q_l\) for each particle (filled with NaN for particles with no neighbors).

  • ave_Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{Q_l}\) for each particle (filled with NaN for particles with no neighbors).

  • norm_Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(Q_l\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

  • ave_norm_Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{Q_l}\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

compute

Compute the order parameter.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAve

Compute the order parameter over two nearest neighbor shells.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAveNorm

Compute the order parameter over two nearest neighbor shells normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeNorm

Compute the order parameter normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

setBox

Reset the simulation box.

Parameters

box (freud.box.Box) – Simulation box.

class freud.order.LocalQlNear(box, rmax, l, kn)

A variant of the LocalQl class that performs its average over nearest neighbor particles as determined by an instance of freud.locality.NeighborList. The number of included neighbors is determined by the kn parameter to the constructor.

Module author: Xiyu Du <xiyudu@umich.edu>

Module author: Vyas Ramasubramani <vramasub@umich.edu>

Parameters
  • box (freud.box.Box) – Simulation box.

  • rmax (float) – Cutoff radius for the local order parameter. Values near the first minimum of the RDF are recommended.

  • l (unsigned int) – Spherical harmonic quantum number l. Must be a positive number.

  • kn (unsigned int) – Number of nearest neighbors. must be a positive integer.

Variables
  • box (freud.box.Box) – Box used in the calculation.

  • num_particles (unsigned int) – Number of particles.

  • Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(Q_l\) for each particle (filled with NaN for particles with no neighbors).

  • ave_Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{Q_l}\) for each particle (filled with NaN for particles with no neighbors).

  • norm_Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(Q_l\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

  • ave_norm_Ql (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{Q_l}\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

compute

Compute the order parameter.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAve

Compute the order parameter over two nearest neighbor shells.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAveNorm

Compute the order parameter over two nearest neighbor shells normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeNorm

Compute the order parameter normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

setBox

Reset the simulation box.

Parameters

box (freud.box.Box) – Simulation box.

Steinhardt \(W_l\) Order Parameter

class freud.order.LocalWl(box, rmax, l)

Compute the local Steinhardt [Steinhardt1983] rotationally invariant \(W_l\) order parameter for a set of points.

Implements the local rotationally invariant \(W_l\) order parameter described by Steinhardt. For a particle i, we calculate the average \(W_l\) by summing the spherical harmonics between particle \(i\) and its neighbors \(j\) in a local region: \(\overline{Q}_{lm}(i) = \frac{1}{N_b} \displaystyle\sum_{j=1}^{N_b} Y_{lm}(\theta(\vec{r}_{ij}), \phi(\vec{r}_{ij}))\). The particles included in the sum are determined by the rmax argument to the constructor.

The \(W_l\) is then defined as a weighted average over the \(\overline{Q}_{lm}(i)\) values using Wigner 3j symbols (Clebsch-Gordan coefficients). The resulting combination is rotationally (i.e. frame) invariant.

The computeAve() method provides access to a variant of this parameter that performs a average over the first and second shell combined [Lechner2008]. To compute this parameter, we perform a second averaging over the first neighbor shell of the particle to implicitly include information about the second neighbor shell. This averaging is performed by replacing the value \(\overline{Q}_{lm}(i)\) in the original definition by the average value of \(\overline{Q}_{lm}(k)\) over all the \(k\) neighbors of particle \(i\) as well as itself.

The computeNorm() and computeAveNorm() methods provide normalized versions of compute() and computeAve(), where the normalization is performed by dividing by the average \(Q_{lm}\) values over all particles.

Module author: Xiyu Du <xiyudu@umich.edu>

Module author: Vyas Ramasubramani <vramasub@umich.edu>

Parameters
  • box (freud.box.Box) – Simulation box.

  • rmax (float) – Cutoff radius for the local order parameter. Values near the first minimum of the RDF are recommended.

  • l (unsigned int) – Spherical harmonic quantum number l. Must be a positive integer.

  • rmin (float) – Lower bound for computing the local order parameter. Allows looking at, for instance, only the second shell, or some other arbitrary RDF region (Default value = 0).

Variables
  • box (freud.box.Box) – Box used in the calculation.

  • num_particles (unsigned int) – Number of particles.

  • Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(W_l\) for each particle (filled with NaN for particles with no neighbors).

  • ave_Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{W}_l\) for each particle (filled with NaN for particles with no neighbors).

  • norm_Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(W_l\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

  • ave_norm_Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{W}_l\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

compute

Compute the order parameter.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAve

Compute the order parameter over two nearest neighbor shells.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAveNorm

Compute the order parameter over two nearest neighbor shells normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeNorm

Compute the order parameter normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

setBox

Reset the simulation box.

Parameters

box (freud.box.Box) – Simulation box.

class freud.order.LocalWlNear(box, rmax, l, kn)

A variant of the LocalWl class that performs its average over nearest neighbor particles as determined by an instance of freud.locality.NeighborList. The number of included neighbors is determined by the kn parameter to the constructor.

Module author: Xiyu Du <xiyudu@umich.edu>

Module author: Vyas Ramasubramani <vramasub@umich.edu>

Parameters
  • box (freud.box.Box) – Simulation box.

  • rmax (float) – Cutoff radius for the local order parameter. Values near the first minimum of the RDF are recommended.

  • l (unsigned int) – Spherical harmonic quantum number l. Must be a positive number

  • kn (unsigned int) – Number of nearest neighbors. Must be a positive number.

Variables
  • box (freud.box.Box) – Box used in the calculation.

  • num_particles (unsigned int) – Number of particles.

  • Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(W_l\) for each particle (filled with NaN for particles with no neighbors).

  • ave_Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{W}_l\) for each particle (filled with NaN for particles with no neighbors).

  • norm_Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(W_l\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

  • ave_norm_Wl (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(\bar{W}_l\) for each particle normalized by the value over all particles (filled with NaN for particles with no neighbors).

compute

Compute the order parameter.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAve

Compute the order parameter over two nearest neighbor shells.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeAveNorm

Compute the order parameter over two nearest neighbor shells normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeNorm

Compute the order parameter normalized by the average spherical harmonic value over all the particles.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

setBox

Reset the simulation box.

Parameters

box (freud.box.Box) – Simulation box.

Solid-Liquid Order Parameter

class freud.order.SolLiq(box, rmax, Qthreshold, Sthreshold, l)

Uses dot products of \(Q_{lm}\) between particles for clustering.

Module author: Richmond Newman <newmanrs@umich.edu>

Parameters
  • box (freud.box.Box) – Simulation box.

  • rmax (float) – Cutoff radius for the local order parameter. Values near first minimum of the RDF are recommended.

  • Qthreshold (float) – Value of dot product threshold when evaluating \(Q_{lm}^*(i) Q_{lm}(j)\) to determine if a neighbor pair is a solid-like bond. (For \(l=6\), 0.7 generally good for FCC or BCC structures).

  • Sthreshold (unsigned int) – Minimum required number of adjacent solid-link bonds for a particle to be considered solid-like for clustering. (For \(l=6\), 6-8 is generally good for FCC or BCC structures).

  • l (unsigned int) – Choose spherical harmonic \(Q_l\). Must be positive and even.

Variables
  • box (freud.box.Box) – Box used in the calculation.

  • largest_cluster_size (unsigned int) – The largest cluster size. Must call a compute method first.

  • cluster_sizes (unsigned int) – The sizes of all clusters.

  • largest_cluster_size – The largest cluster size. Must call a compute method first.

  • Ql_mi (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(Q_{lmi}\) for each particle.

  • clusters (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed set of solid-like cluster indices for each particle.

  • num_connections (\(\left(N_{particles}\right)\) numpy.ndarray) – The number of connections per particle.

  • Ql_dot_ij (\(\left(N_{particles}\right)\) numpy.ndarray) – Reference to the qldot_ij values.

  • num_particles (unsigned int) – Number of particles.

compute

Compute the solid-liquid order parameter.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeSolLiqNoNorm

Compute the solid-liquid order parameter without normalizing the dot product.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

computeSolLiqVariant

Compute a variant of the solid-liquid order parameter.

This variant method places a minimum threshold on the number of solid-like bonds a particle must have to be considered solid-like for clustering purposes.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Points to calculate the order parameter.

  • nlist (freud.locality.NeighborList, optional) – Neighborlist to use to find bonds (Default value = None).

class freud.order.SolLiqNear(box, rmax, Qthreshold, Sthreshold, l, kn)

A variant of the SolLiq class that performs its average over nearest neighbor particles as determined by an instance of freud.locality.NeighborList. The number of included neighbors is determined by the kn parameter to the constructor.

Module author: Richmond Newman <newmanrs@umich.edu>

Parameters
  • box (freud.box.Box) – Simulation box.

  • rmax (float) – Cutoff radius for the local order parameter. Values near the first minimum of the RDF are recommended.

  • Qthreshold (float) – Value of dot product threshold when evaluating \(Q_{lm}^*(i) Q_{lm}(j)\) to determine if a neighbor pair is a solid-like bond. (For \(l=6\), 0.7 generally good for FCC or BCC structures).

  • Sthreshold (unsigned int) – Minimum required number of adjacent solid-link bonds for a particle to be considered solid-like for clustering. (For \(l=6\), 6-8 is generally good for FCC or BCC structures).

  • l (unsigned int) – Choose spherical harmonic \(Q_l\). Must be positive and even.

  • kn (unsigned int) – Number of nearest neighbors. Must be a positive number.

Variables
  • box (freud.box.Box) – Box used in the calculation.

  • largest_cluster_size (unsigned int) – The largest cluster size. Must call a compute method first.

  • cluster_sizes (unsigned int) – The sizes of all clusters.

  • largest_cluster_size – The largest cluster size. Must call a compute method first.

  • Ql_mi (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed \(Q_{lmi}\) for each particle.

  • clusters (\(\left(N_{particles}\right)\) numpy.ndarray) – The last computed set of solid-like cluster indices for each particle.

  • num_connections (\(\left(N_{particles}\right)\) numpy.ndarray) – The number of connections per particle.

  • Ql_dot_ij (\(\left(N_{particles}\right)\) numpy.ndarray) – Reference to the qldot_ij values.

  • num_particles (unsigned int) – Number of particles.

compute

Compute the local rotationally invariant \(Q_l\) order parameter.

Parameters
computeSolLiqNoNorm

Compute the local rotationally invariant \(Q_l\) order parameter.

Parameters
computeSolLiqVariant

Compute the local rotationally invariant \(Q_l\) order parameter.

Parameters

Rotational Autocorrelation

class freud.order.RotationalAutocorrelation(l)

Calculates a measure of total rotational autocorrelation based on hyperspherical harmonics as laid out in “Design rules for engineering colloidal plastic crystals of hard polyhedra - phase behavior and directional entropic forces” by Karas et al. (currently in preparation). The output is not a correlation function, but rather a scalar value that measures total system orientational correlation with an initial state. As such, the output can be treated as an order parameter measuring degrees of rotational (de)correlation. For analysis of a trajectory, the compute call needs to be done at each trajectory frame.

Module author: Andrew Karas <askaras@umich.edu>

Module author: Vyas Ramasubramani <vramasub@umich.edu>

New in version 1.0.

Parameters

l (int) – Order of the hyperspherical harmonic. Must be a positive, even integer.

Variables
  • num_orientations (unsigned int) – The number of orientations used in computing the last set.

  • azimuthal (int) – The azimuthal quantum number, which defines the order of the hyperspherical harmonic. Must be a positive, even integer.

  • ra_array ((\(N_{orientations}\)) numpy.ndarray) – The per-orientation array of rotational autocorrelation values calculated by the last call to compute.

  • autocorrelation (float) – The autocorrelation computed in the last call to compute.

compute

Calculates the rotational autocorrelation function for a single frame.

Parameters
  • ref_ors ((\(N_{orientations}\), 4) numpy.ndarray) – Reference orientations for the initial frame.

  • ors ((\(N_{orientations}\), 4) numpy.ndarray) – Orientations for the frame of interest.