Density Module¶
The density module contains functions which deal with the density of the system
Correlation Functions¶
-
class
freud.density.
FloatCF
(rmax, dr)¶ Computes the pairwise correlation function \(\left< p*q \right> \left( r \right)\) between two sets of points with associated values p and q.
Two sets of points and two sets of real values associated with those points are given. Computing the correlation function results in an array of the expected (average) product of all values at a given radial distance.
The values of r to compute the correlation function at are controlled by the rmax and dr parameters to the constructor. rmax determines the maximum r at which to compute the correlation function and dr is the step size for each bin.
2D: CorrelationFunction properly handles 2D boxes. As with everything else in freud, 2D points must be passed in as 3 component vectors x,y,0. Failing to set 0 in the third component will lead to undefined behavior.
Self-correlation: It is often the case that we wish to compute the correlation function of a set of points with itself. If given the same arrays for both points and ref_points, we omit accumulating the self-correlation value in the first bin.
Module author: Matthew Spellings <mspells@umich.edu>
Parameters: -
R
¶ return – values of bin centers :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float32
-
RDF
¶ return – expected (average) product of all values at a given radial distance :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float64
-
accumulate
(self, box, ref_points, refValues, points, values, nlist=None)¶ Calculates the correlation function and adds to the current histogram.
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – reference points to calculate the local density - refValues (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float64
) – values to use in computation - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – points to calculate the local density - values (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float64
) – values to use in computation - nlist (
freud.locality.NeighborList
) –freud.locality.NeighborList
object to use to find bonds
- box (
-
box
¶ Get the box used in the calculation
Returns: freud Box Return type: freud.box.Box
-
compute
(self, box, ref_points, refValues, points, values, nlist=None)¶ Calculates the correlation function for the given points. Will overwrite the current histogram.
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – reference points to calculate the local density - refValues (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float64
) – values to use in computation - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – points to calculate the local density - values (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float64
) – values to use in computation - nlist (
freud.locality.NeighborList
) –freud.locality.NeighborList
object to use to find bonds
- box (
-
counts
¶ return – counts of each histogram bin :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.int32
-
getBox
(self)¶ Get the box used in the calculation
Returns: freud Box Return type: freud.box.Box
-
getCounts
(self)¶ Returns: counts of each histogram bin Return type: numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.int32
-
getR
(self)¶ Returns: values of bin centers Return type: numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float32
-
getRDF
(self)¶ Returns: expected (average) product of all values at a given radial distance Return type: numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float64
-
reduceCorrelationFunction
(self)¶ Reduces the histogram in the values over N processors to a single histogram. This is called automatically by
freud.density.FloatCF.getRDF()
,freud.density.FloatCF.getCounts()
.
-
resetCorrelationFunction
(self)¶ resets the values of the correlation function histogram in memory
-
-
class
freud.density.
ComplexCF
(rmax, dr)¶ Computes the pairwise correlation function \(\left< p*q \right> \left( r \right)\) between two sets of points with associated values \(p\) and \(q\).
Two sets of points and two sets of complex values associated with those points are given. Computing the correlation function results in an array of the expected (average) product of all values at a given radial distance.
The values of \(r\) to compute the correlation function at are controlled by the rmax and dr parameters to the constructor. rmax determines the maximum r at which to compute the correlation function and dr is the step size for each bin.
2D: CorrelationFunction properly handles 2D boxes. As with everything else in freud, 2D points must be passed in as 3 component vectors x,y,0. Failing to set 0 in the third component will lead to undefined behavior.
Self-correlation: It is often the case that we wish to compute the correlation function of a set of points with itself. If given the same arrays for both points and ref_points, we omit accumulating the self-correlation value in the first bin.
Module author: Matthew Spellings <mspells@umich.edu>
Parameters: -
R
¶ return – values of bin centers :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float32
-
RDF
¶ return – expected (average) product of all values at a given radial distance :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float64
-
accumulate
(self, box, ref_points, refValues, points, values, nlist=None)¶ Calculates the correlation function and adds to the current histogram.
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – reference points to calculate the local density - refValues (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.complex128
) – values to use in computation - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – points to calculate the local density - values (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.complex128
) – values to use in computation - nlist (
freud.locality.NeighborList
) –freud.locality.NeighborList
object to use to find bonds
- box (
-
box
¶ Get the box used in the calculation
Returns: freud Box Return type: freud.box.Box
-
compute
(self, box, ref_points, refValues, points, values, nlist=None)¶ Calculates the correlation function for the given points. Will overwrite the current histogram.
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – reference points to calculate the local density - refValues (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.complex128
) – values to use in computation - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – points to calculate the local density - values (
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.complex128
) – values to use in computation - nlist (
freud.locality.NeighborList
) –freud.locality.NeighborList
object to use to find bonds
- box (
-
counts
¶ return – counts of each histogram bin :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.int32
-
getBox
(self)¶ Returns: freud Box Return type: freud.box.Box()
-
getCounts
(self)¶ Returns: counts of each histogram bin Return type: numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.int32
-
getR
(self)¶ Returns: values of bin centers Return type: numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float32
-
getRDF
(self)¶ Returns: expected (average) product of all values at a given radial distance Return type: numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.complex128
-
reduceCorrelationFunction
(self)¶ Reduces the histogram in the values over N processors to a single histogram. This is called automatically by
freud.density.ComplexCF.getRDF()
,freud.density.ComplexCF.getCounts()
.
-
resetCorrelationFunction
(self)¶ 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 the grid based upon the distance of the grid cell from the center of the Gaussian. The dimensions of the image (grid) are set in the constructor.
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
Constructor Calls:
Initialize with all dimensions identical:
freud.density.GaussianDensity(width, r_cut, dr)
Initialize with each dimension specified:
freud.density.GaussianDensity(width_x, width_y, width_z, r_cut, dr)
-
box
¶ Get the box used in the calculation
Returns: freud Box Return type: freud.box.Box
-
compute
(self, box, points)¶ Calculates the gaussian blur for the specified points. Does not accumulate (will overwrite current image).
Parameters: - box (
freud.box.Box
) – simulation box - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – points to calculate the local density
- box (
-
gaussian_density
¶ return – Image (grid) with values of gaussian :rtype:
numpy.ndarray
, shape=(\(w_x\), \(w_y\), \(w_z\)), dtype=numpy.float32
-
getBox
(self)¶ Returns: freud Box Return type: freud.box.Box
-
getGaussianDensity
(self)¶ Returns: Image (grid) with values of gaussian Return type: numpy.ndarray
, shape=(\(w_x\), \(w_y\), \(w_z\)), dtype=numpy.float32
-
resetDensity
(self)¶ resets the values of GaussianDensity in memory
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 to calculate the local density. volume is the volume of a single particle. diameter is the diameter of the circumsphere of an individual particle.
2D: RDF properly handles 2D boxes. Requires the points to be passed in [x, y, 0]. Failing to z=0 will lead to undefined behavior.
Module author: Joshua Anderson <joaander@umich.edu>
Parameters: -
box
¶ Get the box used in the calculation
Returns: freud Box Return type: freud.box.Box
-
compute
(self, box, ref_points, points=None, nlist=None)¶ Calculates the local density for the specified points. Does not accumulate (will overwrite current data).
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – reference points to calculate the local density - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – (optional) points to calculate the local density - nlist (
freud.locality.NeighborList
) –freud.locality.NeighborList
object to use to find bonds
- box (
-
density
¶ return – Density array for each particle :rtype:
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float32
-
getBox
(self)¶ Returns: freud Box Return type: freud.box.Box
-
getDensity
(self)¶ Returns: Density array for each particle Return type: numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float32
-
getNumNeighbors
(self)¶ Returns: Number of neighbors for each particle Return type: numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float32
-
num_neighbors
¶ return – Number of neighbors for each particle :rtype:
numpy.ndarray
, shape=(\(N_{particles}\)), dtype=numpy.float32
-
Radial Distribution Function¶
-
class
freud.density.
RDF
(rmax, dr)¶ 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 rmax, dr in the constructor. rmax sets the maximum distance at which to calculate the \(g \left( r \right)\) while dr determines the step size for each bin.
Module author: Eric Harper <harperic@umich.edu>
Note
2D: RDF properly handles 2D boxes. Requires the points to be passed in [x, y, 0]. Failing to z=0 will lead to undefined behavior.
Parameters: -
R
¶ return – values of bin centers :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float32
-
RDF
¶ return – expected (average) product of all values at a given radial distance :rtype:
numpy.ndarray
, shape=(\(N_{bins}\)), dtype=numpy.float64
-
accumulate
(self, box, ref_points, points, nlist=None)¶ Calculates the rdf and adds to the current rdf histogram.
Parameters: - box (
freud.box.Box
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – reference points to calculate the local density - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – points to calculate the local density - nlist (
freud.locality.NeighborList
) –freud.locality.NeighborList
object to use to find bonds
- box (
-
box
¶ Get the box used in the calculation
Returns: freud Box Return type: freud.box.Box
-
compute
(self, box, ref_points, points, nlist=None)¶ Calculates the rdf for the specified points. Will overwrite the current histogram.
Parameters: - box (
freud.box.Box()
) – simulation box - ref_points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – reference points to calculate the local density - points (
numpy.ndarray
, shape=(\(N_{particles}\), 3), dtype=numpy.float32
) – points to calculate the local density - nlist (
freud.locality.NeighborList
) –freud.locality.NeighborList
object to use to find bonds
- box (
-
getBox
(self)¶ Returns: freud Box Return type: freud.box.Box
-
getNr
(self)¶ Returns: histogram of cumulative rdf values Return type: numpy.ndarray
, shape=(\(N_{bins}\), 3), dtype=numpy.float32
-
getR
(self)¶ Returns: values of the histogram bin centers Return type: numpy.ndarray
, shape=(\(N_{bins}\), 3), dtype=numpy.float32
-
getRDF
(self)¶ Returns: histogram of rdf values Return type: numpy.ndarray
, shape=(\(N_{bins}\), 3), dtype=numpy.float32
-
n_r
¶ return – histogram of cumulative rdf values :rtype:
numpy.ndarray
, shape=(\(N_{bins}\), 3), dtype=numpy.float32
-
reduceRDF
(self)¶ Reduces the histogram in the values over N processors to a single histogram. This is called automatically by
freud.density.RDF.getRDF()
,freud.density.RDF.getNr()
.
-
resetRDF
(self)¶ resets the values of RDF in memory
-