Density Module¶
Overview
Computes the real pairwise correlation function. |
|
Computes the complex pairwise correlation function. |
|
Computes the density of a system on a grid. |
|
Computes the local density around a particle. |
|
Computes RDF for supplied data. |
Details
The freud.density
module contains various classes relating to the
density of the system. These functions allow evaluation of particle
distributions with respect to other particles.
Correlation Functions¶
-
class
freud.density.
FloatCF
(rmax, dr)¶ Computes the real pairwise correlation function.
The correlation function is given by \(C(r) = \left\langle s_1(0) \cdot s_2(r) \right\rangle\) between two sets of points \(p_1\) (
ref_points
) and \(p_2\) (points
) with associated values \(s_1\) (ref_values
) and \(s_2\) (values
). Computing the correlation function results in an array of the expected (average) product of all values at a given radial distance \(r\).The values of \(r\) where the correlation function is computed are controlled by the
rmax
anddr
parameters to the constructor.rmax
determines the maximum distance at which to compute the correlation function anddr
is the step size for each bin.Note
2D:
freud.density.FloatCF
properly handles 2D boxes. The points must be passed in as[x, y, 0]
. Failing to set z=0 will lead to undefined behavior.Note
Self-correlation: It is often the case that we wish to compute the correlation function of a set of points with itself. If
points
is the same asref_points
, not provided, orNone
, we omit accumulating the self-correlation value in the first bin.Module author: Matthew Spellings <mspells@umich.edu>
- Parameters
- Variables
RDF ((\(N_{bins}\))
numpy.ndarray
) – Expected (average) product of all values whose radial distance falls within a given distance bin.box (
freud.box.Box
) – The box used in the calculation.counts ((\(N_{bins}\))
numpy.ndarray
) – The number of points in each histogram bin.R ((\(N_{bins}\))
numpy.ndarray
) – The centers of each bin.
-
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 the correlation function.ref_values ((\(N_{ref\_points}\))
numpy.ndarray
) – Real values used to calculate the correlation function.points ((\(N_{points}\), 3)
numpy.ndarray
, optional) – Points used to calculate the correlation function. Usesref_points
if not provided orNone
.values ((\(N_{points}\))
numpy.ndarray
, optional) – Real values used to calculate the correlation function. Usesref_values
if not provided orNone
.nlist (
freud.locality.NeighborList
, optional) – NeighborList to use to find bonds (Default value =None
).
-
compute
¶ Calculates the correlation function for the given points. Will overwrite the current histogram.
- Parameters
box (
freud.box.Box
) – Simulation box.ref_points ((\(N_{ref\_points}\), 3)
numpy.ndarray
) – Reference points used to calculate the correlation function.ref_values ((\(N_{ref\_points}\))
numpy.ndarray
) – Real values used to calculate the correlation function.points ((\(N_{points}\), 3)
numpy.ndarray
, optional) – Points used to calculate the correlation function. Usesref_points
if not provided orNone
.values ((\(N_{points}\))
numpy.ndarray
, optional) – Real values used to calculate the correlation function. Usesref_values
if not provided orNone
.nlist (
freud.locality.NeighborList
, optional) – NeighborList to use to find bonds (Default value =None
).
-
plot
¶ Plot correlation function.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. (Default value =None
)- Returns
Axis with the plot.
- Return type
-
reset
¶ Resets the values of the correlation function histogram in memory.
-
class
freud.density.
ComplexCF
(rmax, dr)¶ Computes the complex pairwise correlation function.
The correlation function is given by \(C(r) = \left\langle s_1(0) \cdot s_2(r) \right\rangle\) between two sets of points \(p_1\) (
ref_points
) and \(p_2\) (points
) with associated values \(s_1\) (ref_values
) and \(s_2\) (values
). Computing the correlation function results in an array of the expected (average) product of all values at a given radial distance \(r\).The values of \(r\) where the correlation function is computed are controlled by the
rmax
anddr
parameters to the constructor.rmax
determines the maximum distance at which to compute the correlation function anddr
is the step size for each bin.Note
2D:
freud.density.ComplexCF
properly handles 2D boxes. The points must be passed in as[x, y, 0]
. Failing to set z=0 will lead to undefined behavior.Note
Self-correlation: It is often the case that we wish to compute the correlation function of a set of points with itself. If
points
is the same asref_points
, not provided, orNone
, we omit accumulating the self-correlation value in the first bin.Module author: Matthew Spellings <mspells@umich.edu>
- Parameters
- Variables
RDF ((\(N_{bins}\))
numpy.ndarray
) – Expected (average) product of all values at a given radial distance.box (
freud.box.Box
) – Box used in the calculation.counts ((\(N_{bins}\))
numpy.ndarray
) – The number of points in each histogram bin.R ((\(N_{bins}\))
numpy.ndarray
) – The centers of each bin.
-
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 the correlation function.ref_values ((\(N_{ref\_points}\))
numpy.ndarray
) – Complex values used to calculate the correlation function.points ((\(N_{points}\), 3)
numpy.ndarray
, optional) – Points used to calculate the correlation function. Usesref_points
if not provided orNone
.values ((\(N_{points}\))
numpy.ndarray
, optional) – Complex values used to calculate the correlation function. Usesref_values
if not provided orNone
.nlist (
freud.locality.NeighborList
, optional) – NeighborList to use to find bonds (Default value =None
).
-
compute
¶ Calculates the correlation function for the given points. Will overwrite the current histogram.
- Parameters
box (
freud.box.Box
) – Simulation box.ref_points ((\(N_{ref\_points}\), 3)
numpy.ndarray
) – Reference points used to calculate the correlation function.ref_values ((\(N_{ref\_points}\))
numpy.ndarray
) – Complex values used to calculate the correlation function.points ((\(N_{points}\), 3)
numpy.ndarray
, optional) – Points used to calculate the correlation function. Usesref_points
if not provided orNone
.values ((\(N_{points}\))
numpy.ndarray
, optional) – Complex values used to calculate the correlation function. Usesref_values
if not provided orNone
.nlist (
freud.locality.NeighborList
, optional) – NeighborList to use to find bonds (Default value =None
).
-
plot
¶ Plot complex correlation function.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. (Default value =None
)- Returns
Axis with the plot.
- Return type
-
reset
¶ Resets the values of the correlation function histogram in memory.
Gaussian Density¶
-
class
freud.density.
GaussianDensity
(*args)¶ Computes the density of a system on a grid.
Replaces particle positions with a Gaussian blur and calculates the contribution from each to the proscribed grid based upon the distance of the grid cell from the center of the Gaussian. The resulting data is a regular grid of particle densities that can be used in standard algorithms requiring evenly spaced point, such as Fast Fourier Transforms. The dimensions of the image (grid) are set in the constructor, and can either be set equally for all dimensions or for each dimension independently.
Constructor Calls:
Initialize with all dimensions identical:
freud.density.GaussianDensity(width, r_cut, sigma)
Initialize with each dimension specified:
freud.density.GaussianDensity(width_x, width_y, width_z, r_cut, sigma)
Module author: Joshua Anderson <joaander@umich.edu>
- Parameters
width (unsigned int) – Number of pixels to make the image.
width_x (unsigned int) – Number of pixels to make the image in x.
width_y (unsigned int) – Number of pixels to make the image in y.
width_z (unsigned int) – Number of pixels to make the image in z.
r_cut (float) – Distance over which to blur.
sigma (float) – Sigma parameter for Gaussian.
- Variables
box (
freud.box.Box
) – Box used in the calculation.gaussian_density ((\(w_x\), \(w_y\), \(w_z\))
numpy.ndarray
) – The image grid with the Gaussian density.
-
compute
¶ Calculates the Gaussian blur for the specified points. Does not accumulate (will overwrite current image).
- Parameters
box (
freud.box.Box
) – Simulation box.points ((\(N_{points}\), 3)
numpy.ndarray
) – Points to calculate the local density.
-
plot
¶ Plot Gaussian Density.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. (Default value =None
)- Returns
Axis with the plot.
- Return type
Local Density¶
-
class
freud.density.
LocalDensity
(r_cut, volume, diameter)¶ Computes the local density around a particle.
The density of the local environment is computed and averaged for a given set of reference points in a sea of data points. Providing the same points calculates them against themselves. Computing the local density results in an array listing the value of the local density around each reference point. Also available is the number of neighbors for each reference point, giving the user the ability to count the number of particles in that region.
The values to compute the local density are set in the constructor.
r_cut
sets the maximum distance at which data points are included relative to a given reference point.volume
is the volume of a single data points, anddiameter
is the diameter of the circumsphere of an individual data point. Note that the volume and diameter do not affect the reference point; whether or not data points are counted as neighbors of a given reference point is entirely determined by the distance between reference point and data point center relative tor_cut
and thediameter
of the data point.In order to provide sufficiently smooth data, data points can be fractionally counted towards the density. Rather than perform compute-intensive area (volume) overlap calculations to determine the exact amount of overlap area (volume), the LocalDensity class performs a simple linear interpolation relative to the centers of the data points. Specifically, a point is counted as one neighbor of a given reference point if it is entirely contained within the
r_cut
, half of a neighbor if the distance to its center is exactlyr_cut
, and zero if its center is a distance greater than or equal tor_cut + diameter
from the reference point’s center. Graphically, this looks like:Note
2D:
freud.density.LocalDensity
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: Joshua Anderson <joaander@umich.edu>
- Parameters
- Variables
box (
freud.box.Box
) – Box used in the calculation.density ((\(N_{ref\_points}\))
numpy.ndarray
) – Density of points per ref_point.num_neighbors ((\(N_{ref\_points}\))
numpy.ndarray
) – Number of neighbor points for each ref_point.
-
compute
¶ Calculates the local density for the specified points. Does not accumulate (will overwrite current data).
- Parameters
box (
freud.box.Box
) – Simulation box.ref_points ((\(N_{ref\_points}\), 3)
numpy.ndarray
) – Reference points to calculate the local density.points ((\(N_{points}\), 3)
numpy.ndarray
, optional) – Points to calculate the local density. Usesref_points
if not provided orNone
.nlist (
freud.locality.NeighborList
, optional) – NeighborList to use to find bonds (Default value =None
).
Radial Distribution Function¶
-
class
freud.density.
RDF
(rmax, dr, rmin=0)¶ Computes RDF for supplied data.
The RDF (\(g \left( r \right)\)) is computed and averaged for a given set of reference points in a sea of data points. Providing the same points calculates them against themselves. Computing the RDF results in an RDF array listing the value of the RDF at each given \(r\), listed in the
R
array.The values of \(r\) to compute the RDF are set by the values of
rmin
,rmax
,dr
in the constructor.rmax
sets the maximum distance at which to calculate the \(g \left( r \right)\),rmin
sets the minimum distance at which to calculate the \(g \left( r \right)\), anddr
determines the step size for each bin.Module author: Eric Harper <harperic@umich.edu>
Note
2D:
freud.density.RDF
properly handles 2D boxes. The points must be passed in as[x, y, 0]
. Failing to set z=0 will lead to undefined behavior.- Parameters
- Variables
box (
freud.box.Box
) – Box used in the calculation.RDF ((\(N_{bins}\),)
numpy.ndarray
) – Histogram of RDF values.R ((\(N_{bins}\))
numpy.ndarray
) – The centers of each bin.n_r ((\(N_{bins}\),)
numpy.ndarray
) – Histogram of cumulative RDF values (i.e. the integrated RDF).
Changed in version 0.7.0: Added optional rmin argument.
-
accumulate
¶ Calculates the RDF and adds to the current RDF histogram.
- Parameters
box (
freud.box.Box
) – Simulation box.ref_points ((\(N_{ref\_points}\), 3)
numpy.ndarray
) – Reference points used to calculate the RDF.points ((\(N_{points}\), 3)
numpy.ndarray
, optional) – Points used to calculate the RDF. Usesref_points
if not provided orNone
.nlist (
freud.locality.NeighborList
, optional) – NeighborList to use to find bonds (Default value =None
).
-
compute
¶ Calculates the RDF for the specified points. Will overwrite the current histogram.
- Parameters
box (
freud.box.Box
) – Simulation box.ref_points ((\(N_{ref\_points}\), 3)
numpy.ndarray
) – Reference points used to calculate the RDF.points ((\(N_{points}\), 3)
numpy.ndarray
, optional) – Points used to calculate the RDF. Usesref_points
if not provided orNone
.nlist (
freud.locality.NeighborList
) – NeighborList to use to find bonds (Default value =None
).
-
plot
¶ Plot radial distribution function.
- Parameters
ax (
matplotlib.axes.Axes
) – Axis to plot on. IfNone
, make a new figure and axis. (Default value =None
)- Returns
Axis with the plot.
- Return type
-
reset
¶ Resets the values of RDF in memory.