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.EnvironmentCluster

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

freud.environment.EnvironmentMotifMatch

Find matches between local arrangements of a set of points and a provided motif.

freud.environment.AngularSeparationGlobal

Calculates the minimum angles of separation between orientations and global orientations.

freud.environment.AngularSeparationNeighbor

Calculates the minimum angles of separation between orientations and query orientations.

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.

class freud.environment.AngularSeparationGlobal

Bases: freud.util._Compute

Calculates the minimum angles of separation between orientations and global orientations.

property angles

\(\left(N_{orientations}, N_{global\_orientations}\right)\) numpy.ndarray: The global angles in radians.

compute

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

Parameters
  • global_orientations ((\(N_{global}\), 4) numpy.ndarray) – Set of global orientations to calculate the order parameter.

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

  • equiv_orientations ((\(N_{equiv}\), 4) numpy.ndarray, optional) – The set of all equivalent quaternions that map the particle to itself (the elements of its rotational symmetry group). Important: equiv_orientations must include both \(q\) and \(-q\), for all included quaternions. Note that this calculation assumes that all points in the system share the same set of equivalent orientations. (Default value = [[1, 0, 0, 0]])

class freud.environment.AngularSeparationNeighbor

Bases: freud.locality._PairCompute

Calculates the minimum angles of separation between orientations and query orientations.

property angles

The neighbor angles in radians. The angles are stored in the order of the neighborlist object.

Type

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

compute

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

Parameters
  • system – Any object that is a valid argument to freud.locality.NeighborQuery.from_system.

  • orientations ((\(N_{points}\), 4) numpy.ndarray) – Orientations associated with system points that are used to calculate bonds.

  • query_points ((\(N_{query\_points}\), 3) numpy.ndarray, optional) – Query points used to calculate the correlation function. Uses the system’s points if None (Default value = None).

  • query_orientations ((\(N_{query\_points}\), 4) numpy.ndarray, optional) – Query orientations used to calculate bonds. Uses orientations if None. (Default value = None).

  • equiv_orientations ((\(N_{equiv}\), 4) numpy.ndarray, optional) – The set of all equivalent quaternions that map the particle to itself (the elements of its rotational symmetry group). Important: equiv_orientations must include both \(q\) and \(-q\), for all included quaternions. Note that this calculation assumes that all points in the system share the same set of equivalent orientations. (Default value = [[1, 0, 0, 0]])

  • neighbors (freud.locality.NeighborList or dict, optional) – Either a NeighborList of neighbor pairs to use in the calculation, or a dictionary of query arguments (Default value: None).

default_query_args

No default query arguments.

property nlist

The neighbor list from the last compute.

Type

freud.locality.NeighborList

class freud.environment.BondOrder

Bases: freud.locality._SpatialHistogram

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() method. 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.

Parameters
  • bins (unsigned int or sequence of length 2) – If an unsigned int, the number of bins in \(\theta\) and \(\phi\). If a sequence of two integers, interpreted as (num_bins_theta, num_bins_phi).

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

bin_centers

The centers of each bin in the histogram (has the same shape as the histogram itself).

Type

numpy.ndarray

property bin_counts

The bin counts in the histogram.

Type

numpy.ndarray

bin_edges

The edges of each bin in the histogram (is one element larger in each dimension than the histogram because each bin has a lower and upper bound).

Type

numpy.ndarray

property bond_order

Bond order.

Type

\(\left(N_{\phi}, N_{\theta} \right)\) numpy.ndarray

bounds

A list of tuples indicating upper and lower bounds of each axis of the histogram.

Type

list (tuple)

property box

Box used in the calculation.

Type

freud.box.Box

compute

Calculates the correlation function and adds to the current histogram.

Parameters
  • system – Any object that is a valid argument to freud.locality.NeighborQuery.from_system.

  • orientations ((\(N_{points}\), 4) numpy.ndarray) – Orientations associated with system points that are used to calculate bonds. Uses identity quaternions if None (Default value = None).

  • query_points ((\(N_{query\_points}\), 3) numpy.ndarray, optional) – Query points used to calculate the correlation function. Uses the system’s points if None (Default value = None).

  • query_orientations ((\(N_{query\_points}\), 4) numpy.ndarray, optional) – Query orientations used to calculate bonds. Uses orientations if None. (Default value = None).

  • neighbors (freud.locality.NeighborList or dict, optional) –

    Either a NeighborList of neighbor pairs to use in the calculation, or a dictionary of query arguments (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).

default_query_args

No default query arguments.

mode

Bond order mode.

Type

str

nbins

The number of bins in each dimension of the histogram

Type

list

class freud.environment.EnvironmentCluster

Bases: freud.environment._MatchEnv

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

property cluster_environments

\(\left(N_{clusters}, N_{neighbors}, 3\right)\) numpy.ndarray): The environments for all clusters.

property cluster_idx

The per-particle index indicating cluster membership.

Type

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

compute

Determine clusters of particles with matching environments.

In general, it is recommended to specify a number of neighbors rather than just a distance cutoff as part of your neighbor querying when performing this computation. Using a distance cutoff alone could easily lead to situations where a point doesn’t match a cluster because a required neighbor is just outside the cutoff.

Parameters
  • system – Any object that is a valid argument to freud.locality.NeighborQuery.from_system.

  • threshold (float) – Maximum magnitude of the vector difference between two vectors, below which they are “matching”. Typically, a good choice is between 10% and 30% of the first well in the radial distribution function (this has distance units).

  • neighbors (freud.locality.NeighborList or dict, optional) –

    Either a NeighborList of neighbor pairs to use in the calculation, or a dictionary of query arguments (Default value: None).

  • env_neighbors (freud.locality.NeighborList or dict, optional) –

    Either a NeighborList of neighbor pairs to use in the calculation, or a dictionary of query arguments (Default value: None). This argument is used to define the neighbors of the environment that motifs are registered against.

  • 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)

  • global_search (bool, optional) – 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. (Default value = False)

default_query_args

No default query arguments.

property num_clusters

The number of clusters.

Type

unsigned int

plot

Plot cluster distribution.

Parameters

ax (matplotlib.axes.Axes, optional) – 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)

property point_environments

\(\left(N_{points}, N_{neighbors}, 3\right)\) numpy.ndarray: All environments for all points.

class freud.environment.EnvironmentMotifMatch

Bases: freud.environment._MatchEnv

Find matches between local arrangements of a set of points and a provided motif.

In general, it is recommended to specify a number of neighbors rather than just a distance cutoff as part of your neighbor querying when performing this computation since it can otherwise be very sensitive. Specifically, it is highly recommended that you choose a number of neighbors query that requests at least as many neighbors as the size of the motif you intend to test against. Otherwise, you will struggle to match the motif. However, this is not currently enforced.

compute

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

Parameters
  • system – Any object that is a valid argument to freud.locality.NeighborQuery.from_system.

  • motif ((\(N_{points}\), 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 “matching”. Typically, a good choice is between 10% and 30% of the first well in the radial distribution function (this has distance units).

  • neighbors (freud.locality.NeighborList or dict, optional) –

    Either a NeighborList of neighbor pairs to use in the calculation, or a dictionary of query arguments (Default value: None).

  • 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).

default_query_args

No default query arguments.

property matches

A boolean array indicating whether each point matches the motif.

Type

\((N_points, )\) numpy.ndarray

property point_environments

\(\left(N_{points}, N_{neighbors}, 3\right)\) numpy.ndarray: All environments for all points.

class freud.environment.LocalBondProjection

Bases: freud.locality._PairCompute

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

compute

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

Parameters
  • system – Any object that is a valid argument to freud.locality.NeighborQuery.from_system.

  • orientations ((\(N_{points}\), 4) numpy.ndarray) – Orientations associated with system points that are used to calculate bonds.

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

  • query_points ((\(N_{query\_points}\), 3) numpy.ndarray, optional) – Query points used to calculate the correlation function. Uses the system’s points if None (Default value = None). (Default value = None).

  • equiv_orientations ((\(N_{equiv}\), 4) numpy.ndarray, optional) – The set of all equivalent quaternions that map the particle to itself (the elements of its rotational symmetry group). Important: equiv_orientations must include both \(q\) and \(-q\), for all included quaternions. Note that this calculation assumes that all points in the system share the same set of equivalent orientations. (Default value = [[1, 0, 0, 0]])

  • neighbors (freud.locality.NeighborList or dict, optional) –

    Either a NeighborList of neighbor pairs to use in the calculation, or a dictionary of query arguments (Default value: None).

default_query_args

No default query arguments.

property nlist

The neighbor list from the last compute.

Type

freud.locality.NeighborList

property normed_projections

\(\left(N_{bonds}, N_{projection\_vecs} \right)\) numpy.ndarray: The projection of each bond between query particles and their neighbors onto each of the projection vectors, normalized by the length of the bond.

property projections

\(\left(N_{bonds}, N_{projection\_vecs} \right)\) numpy.ndarray: The projection of each bond between query particles and their neighbors onto each of the projection vectors.

class freud.environment.LocalDescriptors

Bases: freud.locality._PairCompute

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.

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

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

  • 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').

compute

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

Parameters
default_query_args

No default query arguments.

l_max

The maximum spherical harmonic \(l\) calculated for.

Type

unsigned int

mode

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.

Type

str

negative_m

True if we also calculated \(Y_{lm}\) for negative \(m\).

Type

bool

property nlist

The neighbor list from the last compute.

Type

freud.locality.NeighborList

property num_sphs

The last number of spherical harmonics computed. This is equal to the number of bonds in the last computation, which is at most the number of points multiplied by the lower of the num_neighbors arguments passed to the last compute call or the constructor (it may be less if there are not enough neighbors for every particle).

Type

unsigned int

property sph

\(\left(N_{bonds}, \text{SphWidth} \right)\) numpy.ndarray: The last computed spherical harmonic array.