Environment Module

Overview

freud.environment.BondOrder

Compute the bond orientational order diagram for the system of particles.

freud.environment.LocalDescriptors

Compute a set of descriptors (a numerical “fingerprint”) of a particle’s local environment.

freud.environment.MatchEnv

Clusters particles according to whether their local environments match or not, according to various shape matching metrics.

freud.environment.AngularSeparation

Calculates the minimum angles of separation between particles and references.

freud.environment.LocalBondProjection

Calculates the maximal projection of nearest neighbor bonds for each particle onto some set of reference vectors, defined in the particles’ local reference frame.

Details

The freud.environment module contains functions which characterize the local environments of particles in the system. These methods use the positions and orientations of particles in the local neighborhood of a given particle to characterize the particle environment.

Bond Order

class freud.environment.BondOrder(rmax, k, n, nBinsT, nBinsP)

Compute the bond orientational order diagram for the system of particles.

The bond orientational order diagram (BOOD) is a way of studying the average local environments experienced by particles. In a BOOD, a particle and its nearest neighbors (determined by either a prespecified number of neighbors or simply a cutoff distance) are treated as connected by a bond joining their centers. All of the bonds in the system are then binned by their azimuthal (\(\theta\)) and polar (\(\phi\)) angles to indicate the location of a particle’s neighbors relative to itself. The distance between the particle and its neighbor is only important when determining whether it is counted as a neighbor, but is not part of the BOOD; as such, the BOOD can be viewed as a projection of all bonds onto the unit sphere. The resulting 2D histogram provides insight into how particles are situated relative to one-another in a system.

This class provides access to the classical BOOD as well as a few useful variants. These variants can be accessed via the mode arguments to the compute() or accumulate() methods. Available modes of calculation are:

  • 'bod' (Bond Order Diagram, default): This mode constructs the default BOOD, which is the 2D histogram containing the number of bonds formed through each azimuthal \(\left( \theta \right)\) and polar \(\left( \phi \right)\) angle.

  • 'lbod' (Local Bond Order Diagram): In this mode, a particle’s neighbors are rotated into the local frame of the particle before the BOOD is calculated, i.e. the directions of bonds are determined relative to the orientation of the particle rather than relative to the global reference frame. An example of when this mode would be useful is when a system is composed of multiple grains of the same crystal; the normal BOOD would show twice as many peaks as expected, but using this mode, the bonds would be superimposed.

  • 'obcd' (Orientation Bond Correlation Diagram): This mode aims to quantify the degree of orientational as well as translational ordering. As a first step, the rotation that would align a particle’s neighbor with the particle is calculated. Then, the neighbor is rotated around the central particle by that amount, which actually changes the direction of the bond. One example of how this mode could be useful is in identifying plastic crystals, which exhibit translational but not orientational ordering. Normally, the BOOD for a plastic crystal would exhibit clear structure since there is translational order, but with this mode, the neighbor positions would actually be modified, resulting in an isotropic (disordered) BOOD.

  • 'oocd' (Orientation Orientation Correlation Diagram): This mode is substantially different from the other modes. Rather than compute the histogram of neighbor bonds, this mode instead computes a histogram of the directors of neighboring particles, where the director is defined as the basis vector \(\hat{z}\) rotated by the neighbor’s quaternion. The directors are then rotated into the central particle’s reference frame. This mode provides insight into the local orientational environment of particles, indicating, on average, how a particle’s neighbors are oriented.

Module author: Erin Teich <erteich@umich.edu>

Parameters
  • rmax (float) – Distance over which to calculate.

  • k (unsigned int) – Order parameter i. To be removed.

  • n (unsigned int) – Number of neighbors to find.

  • n_bins_t (unsigned int) – Number of \(\theta\) bins.

  • n_bins_p (unsigned int) – Number of \(\phi\) bins.

Variables
  • bond_order (\(\left(N_{\phi}, N_{\theta} \right)\) numpy.ndarray) – Bond order.

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

  • theta (\(\left(N_{\theta} \right)\) numpy.ndarray) – The values of bin centers for \(\theta\).

  • phi (\(\left(N_{\phi} \right)\) numpy.ndarray) – The values of bin centers for \(\phi\).

  • n_bins_theta (unsigned int) – The number of bins in the \(\theta\) dimension.

  • n_bins_phi (unsigned int) – The number of bins in the \(\phi\) dimension.

accumulate

Calculates the correlation function and adds to the current histogram.

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

  • ref_points ((\(N_{ref\_points}\), 3) numpy.ndarray) – Reference points used to calculate bonds.

  • ref_orientations ((\(N_{ref\_points}\), 4) numpy.ndarray) – Reference orientations used to calculate bonds.

  • points ((\(N_{points}\), 3) numpy.ndarray, optional) – Points used to calculate bonds. Uses ref_points if not provided or None.

  • orientations ((\(N_{points}\), 4) numpy.ndarray) – Orientations used to calculate bonds. Uses ref_orientations if not provided or None.

  • mode (str, optional) – Mode to calculate bond order. Options are 'bod', 'lbod', 'obcd', or 'oocd' (Default value = 'bod').

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

compute

Calculates the bond order histogram. Will overwrite the current histogram.

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

  • ref_points ((\(N_{particles}\), 3) numpy.ndarray) – Reference points used to calculate bonds.

  • ref_orientations ((\(N_{particles}\), 4) numpy.ndarray) – Reference orientations used to calculate bonds.

  • points ((\(N_{particles}\), 3) numpy.ndarray, optional) – Points used to calculate bonds. Uses ref_points if not provided or None.

  • orientations ((\(N_{particles}\), 4) numpy.ndarray, optional) – Orientations used to calculate bonds. Uses ref_orientations if not provided or None.

  • mode (str, optional) – Mode to calculate bond order. Options are 'bod', 'lbod', 'obcd', or 'oocd' (Default value = 'bod').

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

reset

Resets the values of the bond order in memory.

Local Descriptors

class freud.environment.LocalDescriptors(num_neighbors, lmax, rmax, negative_m=True)

Compute a set of descriptors (a numerical “fingerprint”) of a particle’s local environment.

The resulting spherical harmonic array will be a complex-valued array of shape (num_bonds, num_sphs). Spherical harmonic calculation can be restricted to some number of nearest neighbors through the num_neighbors argument; if a particle has more bonds than this number, the last one or more rows of bond spherical harmonics for each particle will not be set.

Module author: Matthew Spellings <mspells@umich.edu>

Parameters
  • num_neighbors (unsigned int) – Maximum number of neighbors to compute descriptors for.

  • lmax (unsigned int) – Maximum spherical harmonic \(l\) to consider.

  • rmax (float) – Initial guess of the maximum radius to looks for neighbors.

  • negative_m (bool) – True if we should also calculate \(Y_{lm}\) for negative \(m\).

Variables
  • sph (\(\left(N_{bonds}, \text{SphWidth} \right)\) numpy.ndarray) – A reference to the last computed spherical harmonic array.

  • num_particles (unsigned int) – The number of points passed to the last call to compute().

  • num_neighbors (unsigned int) – The number of neighbors used by the last call to compute. Bounded from above by the number of reference points multiplied by the lower of the num_neighbors arguments passed to the last compute call or the constructor.

  • l_max (unsigned int) – The maximum spherical harmonic \(l\) to calculate for.

  • r_max (float) – The cutoff radius.

compute

Calculates the local descriptors of bonds from a set of source points to a set of destination points.

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

  • num_neighbors (unsigned int) – Number of nearest neighbors to compute with or to limit to, if the neighbor list is precomputed.

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

  • points ((\(N_{particles}\), 3) numpy.ndarray, optional) – Destination points to calculate the order parameter (Default value = None).

  • orientations ((\(N_{particles}\), 4) numpy.ndarray, optional) – Orientation of each reference point (Default value = None).

  • mode (str, optional) – Orientation mode to use for environments, either 'neighborhood' to use the orientation of the local neighborhood, 'particle_local' to use the given particle orientations, or 'global' to not rotate environments (Default value = 'neighborhood').

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

Match Environments

class freud.environment.MatchEnv(box, rmax, k)

Clusters particles according to whether their local environments match or not, according to various shape matching metrics.

Module author: Erin Teich <erteich@umich.edu>

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

  • rmax (float) – Cutoff radius for cell list and clustering algorithm. Values near the first minimum of the RDF are recommended.

  • k (unsigned int) – Number of nearest neighbors taken to define the local environment of any given particle.

Variables
  • tot_environment (\(\left(N_{particles}, N_{neighbors}, 3\right)\) numpy.ndarray) – All environments for all particles.

  • num_particles (unsigned int) – The number of particles.

  • num_clusters (unsigned int) – The number of clusters.

  • clusters (\(\left(N_{particles}\right)\) numpy.ndarray) – The per-particle index indicating cluster membership.

cluster

Determine clusters of particles with matching environments.

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

  • threshold (float) – Maximum magnitude of the vector difference between two vectors, below which they are “matching.”

  • hard_r (bool) – If True, add all particles that fall within the threshold of m_rmaxsq to the environment.

  • registration (bool) – If True, first use brute force registration to orient one set of environment vectors with respect to the other set such that it minimizes the RMSD between the two sets.

  • global_search (bool) – If True, do an exhaustive search wherein the environments of every single pair of particles in the simulation are compared. If False, only compare the environments of neighboring particles.

  • env_nlist (freud.locality.NeighborList, optional) – NeighborList to use to find the environment of every particle (Default value = None).

  • nlist (freud.locality.NeighborList, optional) – NeighborList to use to find neighbors of every particle, to compare environments (Default value = None).

getEnvironment

Returns the set of vectors defining the environment indexed by i.

Parameters

i (unsigned int) – Environment index.

Returns

The array of vectors.

Return type

\(\left(N_{neighbors}, 3\right)\) numpy.ndarray

isSimilar

Test if the motif provided by refPoints1 is similar to the motif provided by refPoints2.

Parameters
  • refPoints1 ((\(N_{particles}\), 3) numpy.ndarray) – Vectors that make up motif 1.

  • refPoints2 ((\(N_{particles}\), 3) numpy.ndarray) – Vectors that make up motif 2.

  • threshold (float) – Maximum magnitude of the vector difference between two vectors, below which they are considered “matching.”

  • registration (bool, optional) – If True, first use brute force registration to orient one set of environment vectors with respect to the other set such that it minimizes the RMSD between the two sets (Default value = False).

Returns

A doublet that gives the rotated (or not) set of refPoints2, and the mapping between the vectors of refPoints1 and refPoints2 that will make them correspond to each other. Empty if they do not correspond to each other.

Return type

tuple ((\(\left(N_{particles}, 3\right)\) numpy.ndarray), map[int, int])

matchMotif

Determine clusters of particles that match the motif provided by refPoints.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Particle positions.

  • refPoints ((\(N_{particles}\), 3) numpy.ndarray) – Vectors that make up the motif against which we are matching.

  • threshold (float) – Maximum magnitude of the vector difference between two vectors, below which they are considered “matching.”

  • registration (bool, optional) – If True, first use brute force registration to orient one set of environment vectors with respect to the other set such that it minimizes the RMSD between the two sets (Default value = False).

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

minRMSDMotif

Rotate (if registration=True) and permute the environments of all particles to minimize their RMSD with respect to the motif provided by refPoints.

Parameters
  • points ((\(N_{particles}\), 3) numpy.ndarray) – Particle positions.

  • refPoints ((\(N_{particles}\), 3) numpy.ndarray) – Vectors that make up the motif against which we are matching.

  • registration (bool, optional) – If True, first use brute force registration to orient one set of environment vectors with respect to the other set such that it minimizes the RMSD between the two sets (Default value = False).

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

Returns

Vector of minimal RMSD values, one value per particle.

Return type

\(\left(N_{particles}\right)\) numpy.ndarray

minimizeRMSD

Get the somewhat-optimal RMSD between the set of vectors refPoints1 and the set of vectors refPoints2.

Parameters
  • refPoints1 ((\(N_{particles}\), 3) numpy.ndarray) – Vectors that make up motif 1.

  • refPoints2 ((\(N_{particles}\), 3) numpy.ndarray) – Vectors that make up motif 2.

  • registration (bool, optional) – If true, first use brute force registration to orient one set of environment vectors with respect to the other set such that it minimizes the RMSD between the two sets (Default value = False).

Returns

A triplet that gives the associated min_rmsd, rotated (or not) set of refPoints2, and the mapping between the vectors of refPoints1 and refPoints2 that somewhat minimizes the RMSD.

Return type

tuple (float, (\(\left(N_{particles}, 3\right)\) numpy.ndarray), map[int, int])

plot

Plot cluster distribution.

Parameters

ax (matplotlib.axes.Axes) – Axis to plot on. If None, make a new figure and axis. (Default value = None)

Returns

Axis with the plot.

Return type

(matplotlib.axes.Axes)

setBox

Reset the simulation box.

Parameters

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

Angular Separation

class freud.environment.AngularSeparation(rmax, n)

Calculates the minimum angles of separation between particles and references.

Module author: Erin Teich <erteich@umich.edu>

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

Parameters
  • rmax (float) – Cutoff radius for cell list and clustering algorithm. Values near the first minimum of the RDF are recommended.

  • n (int) – The number of neighbors.

Variables
  • nlist (freud.locality.NeighborList) – The neighbor list.

  • n_p (unsigned int) – The number of particles used in computing the last set.

  • n_ref (unsigned int) – The number of reference particles used in computing the neighbor angles.

  • n_global (unsigned int) – The number of global orientations to check against.

  • neighbor_angles (\(\left(N_{bonds}\right)\) numpy.ndarray) – The neighbor angles in radians. This field is only populated after computeNeighbor() is called. The angles are stored in the order of the neighborlist object.

  • global_angles (\(\left(N_{particles}, N_{global} \right)\) numpy.ndarray) – The global angles in radians. This field is only populated after computeGlobal() is called. The angles are stored in the order of the neighborlist object.

computeGlobal

Calculates the minimum angles of separation between global_ors and ors, checking for underlying symmetry as encoded in equiv_quats. The result is stored in the global_angles class attribute.

Parameters
  • global_ors ((\(N_{particles}\), 4) numpy.ndarray) – Reference orientations to calculate the order parameter.

  • ors ((\(N_{particles}\), 4) numpy.ndarray) – Orientations to calculate the order parameter.

  • equiv_quats ((\(N_{particles}\), 4) numpy.ndarray) – The set of all equivalent quaternions that takes the particle as it is defined to some global reference orientation. Important: equiv_quats must include both \(q\) and \(-q\), for all included quaternions.

computeNeighbor

Calculates the minimum angles of separation between ref_ors and ors, checking for underlying symmetry as encoded in equiv_quats. The result is stored in the neighbor_angles class attribute.

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

  • ref_ors ((\(N_{particles}\), 4) numpy.ndarray) – Reference orientations used to calculate the order parameter.

  • ors ((\(N_{particles}\), 4) numpy.ndarray) – Orientations used to calculate the order parameter.

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

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

  • equiv_quats ((\(N_{particles}\), 4) numpy.ndarray, optional) – The set of all equivalent quaternions that takes the particle as it is defined to some global reference orientation. Important: equiv_quats must include both \(q\) and \(-q\), for all included quaternions.

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

Local Bond Projection

class freud.environment.LocalBondProjection(rmax, num_neighbors)

Calculates the maximal projection of nearest neighbor bonds for each particle onto some set of reference vectors, defined in the particles’ local reference frame.

Module author: Erin Teich <erteich@umich.edu>

New in version 0.11.

Parameters
  • rmax (float) – Cutoff radius.

  • num_neighbors (unsigned int) – The number of neighbors.

Variables
  • projections ((\(\left(N_{reference}, N_{neighbors}, N_{projection\_vecs} \right)\) numpy.ndarray) – The projection of each bond between reference particles and their neighbors onto each of the projection vectors.

  • normed_projections ((\(\left(N_{reference}, N_{neighbors}, N_{projection\_vecs} \right)\) numpy.ndarray) – The normalized projection of each bond between reference particles and their neighbors onto each of the projection vectors.

  • num_reference_particles (int) – The number of reference points used in the last calculation.

  • num_particles (int) – The number of points used in the last calculation.

  • num_proj_vectors (int) – The number of projection vectors used in the last calculation.

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

  • nlist (freud.locality.NeighborList) – The neighbor list generated in the last calculation.

compute

Calculates the maximal projections of nearest neighbor bonds (between ref_points and points) onto the set of reference vectors proj_vecs, defined in the local reference frames of the ref_points as defined by the orientations ref_ors. This computation accounts for the underlying symmetries of the reference frame as encoded in equiv_quats.

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

  • proj_vecs ((\(N_{vectors}\), 3) numpy.ndarray) – The set of reference vectors, defined in the reference particles’ reference frame, to calculate maximal local bond projections onto.

  • ref_points ((\(N_{particles}\), 3) numpy.ndarray) – Reference points used in the calculation.

  • ref_ors ((\(N_{particles}\), 4) numpy.ndarray) – Reference orientations used in the calculation.

  • points ((\(N_{particles}\), 3) numpy.ndarray, optional) – Points (neighbors of ref_points) used in the calculation. Uses ref_points if not provided or None.

  • equiv_quats ((\(N_{quats}\), 4) numpy.ndarray, optional) – The set of all equivalent quaternions that takes the particle as it is defined to some global reference orientation. Note that this does not need to include both \(q\) and \(-q\), since \(q\) and \(-q\) effect the same rotation on vectors. Defaults to an identity quaternion.

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