RadialDistributionFunction¶
- class mdcraft.analysis.structure.RadialDistributionFunction(group_a: AtomGroup, group_b: AtomGroup = None, n_bins: int = 201, range_: tuple[float] = (0.0, 15.0), *, exclusion: tuple[int] = None, groupings: str | tuple[str] = 'atoms', drop_axis: int | str = None, norm: str = 'rdf', n_batches: int = None, reduced: bool = False, parallel: bool = False, verbose: bool = True, **kwargs)[source]¶
Bases:
DynamicAnalysisBaseSerial and parallel implementations to calculate the radial distribution function (RDF) \(g_{\alpha\beta}(r)\) between types \(\alpha\) and \(\beta\) and its related properties for two- and three-dimensional systems.
The RDF is given by
\[\begin{split}g_{\alpha\beta}^\mathrm{3D}(r) =\frac{V}{4\pi r^2N_\alpha N_\beta}\sum_{i=1}^{N_\alpha} \sum_{j=1}^{N_\beta}\left\langle\delta \left(|\mathbf{r}_i-\mathbf{r}_j|-r\right)\right\rangle\\ g_{\alpha\beta}^\mathrm{2D}(r) =\frac{A}{2\pi rN_\alpha N_\beta}\sum_{i=1}^{N_\alpha} \sum_{j=1}^{N_\beta}\left\langle\delta \left(|\mathbf{r}_i-\mathbf{r}_j|-r\right)\right\rangle\end{split}\]where \(V\) and \(A\) are the system volume and area, \(N_\alpha\) and \(N_\beta\) are the number of particles of species \(\alpha\) and \(\beta\), and \(\mathbf{r}_i\) and \(\mathbf{r}_j\) are the positions of particles \(i\) and \(j\) belonging to species \(\alpha\) and \(\beta\), respectively. The RDF is normalized such that \(\lim_{r\rightarrow\infty}g_{\alpha\beta}(r)=1\) in a homogeneous system.
A closely related quantity is the single particle density \(n_{\alpha\beta}(r)=\rho_\beta g_{\alpha\beta}(r)\), where \(\rho_\beta\) is the number density of species \(\beta\).
The cumulative RDF is
\[\begin{split}G_{\alpha\beta}^\mathrm{3D}(r) =4\pi\int_0^rR^2g_{\alpha\beta}(R)\,dR\\ G_{\alpha\beta}^\mathrm{2D}(r) =2\pi\int_0^rRg_{\alpha\beta}(R)\,dR\end{split}\]and the average number of \(\beta\) particles found within radius \(r\) is
\[N_{\alpha\beta}(r)=\rho_\beta G_{\alpha\beta}(r)\]The expression above can be used to compute the coordination numbers (number of neighbors in each solvation shell) by setting \(r\) to values where \(g_{\alpha\beta}(r)\) is locally minimized, which signify the solvation shell boundaries.
The RDF can also be used to obtain other relevant structural properties, such as
the potential of mean force
\[w_{\alpha\beta}(r) =-k_\mathrm{B}T\ln{g_{\alpha\beta}(r)}\]where \(k_\mathrm{B}\) is the Boltzmann constant and \(T\) is the system temperature, and
the static or partial structure factor (see
calculate_structure_factor()for the possible definitions).
- Parameters:
- group_aMDAnalysis.AtomGroup
Atom group containing species \(\alpha\).
- group_bMDAnalysis.AtomGroup
Atom group containing species \(\beta\). If not specified, group_a is used for both species.
- n_binsint, default:
201 Number of histogram bins \(N_\mathrm{bins}\).
- range_array-like, default:
(0.0, 15.0) Range of radii values. The upper bound should be less than half the smallest system dimension.
Shape: \((2,)\).
Reference unit: \(\mathrm{Å}\).
- exclusionarray-like, keyword-only, optional
Tiles to exclude from the interparticle distances. The groupings parameter dictates what a tile represents.
Shape: \((2,)\).
Example:
(1, 1)to exclude self-interactions.- groupingsstr or array-like, keyword-only, default:
"atoms" Determines whether the centers of mass are used in lieu of individual atom positions. If groupings is a str, the same value is used for all groups.
Note
If the desired grouping is not
"atoms",the trajectory file should have segments (or chains) containing residues (or molecules) and residues containing atoms, and
residues and segments should be locally unwrapped at the simulation box edges, if not already, using
MDAnalysis.transformations.wrap.unwrap,MDAnalysis.core.groups.AtomGroup.unwrap(), orMDAnalysis.lib.mdamath.make_whole().
Valid values:
"atoms": Atom positions (generally or for coarse-grained simulations)."residues": Residues’ centers of mass (for atomistic simulations)."segments": Segments’ centers of mass (for atomistic polymer simulations).
- drop_axisint or str, keyword-only, optional
Axis in three-dimensional space to ignore in the two-dimensional analysis.
Valid values:
0orxfor the \(x\)-axis,1oryfor the \(y\)-axis, and2orzfor the \(z\)-axis.- normstr, keyword-only, default:
"rdf" Determines how the radial histograms are normalized.
Valid values:
norm="rdf": The radial distribution function \(g_{\alpha\beta}(r)\) is computed.norm="density": The single particle density \(n_{\alpha\beta}(r)\) is computed.norm=None: The raw particle pair count in the radial histogram bins is returned.
- n_batchesint, keyword-only, optional
Number of batches to divide the histogram calculation into. This is useful for large systems that cannot be processed in a single pass.
Note
If you use too few bins and too many batches, the histogram counts may be off by a few due to the floating-point nature of the cutoffs. However, when the RDF is averaged over a long trajectory with many particles, the difference should be negligible.
- reducedbool, keyword-only, default:
False Specifies whether the data is in reduced units.
- parallelbool, keyword-only, default:
False Determines whether the analysis is performed in parallel.
Note
The Dask backend generally provides the best performance.
- verbosebool, keyword-only, default:
True Determines whether detailed progress is shown.
- **kwargs
Additional keyword arguments to pass to
MDAnalysis.analysis.base.AnalysisBase.
- Attributes:
- universeMDAnalysis.Universe
MDAnalysis.core.universe.Universeobject containing all information describing the system.- results.unitsdict
Reference units for the results. For example, to get the reference units for
results.bins, callresults.units["bins"].- results.bin_edgesnumpy.ndarray
Bin edges.
Shape: \((N_\mathrm{bins}+1,)\).
Reference unit: \(\mathrm{Å}\).
- results.binsnumpy.ndarray
Bin centers \(r\).
Shape: \((N_\mathrm{bins},)\).
Reference unit: \(\mathrm{Å}\).
- results.countsnumpy.ndarray
Raw particle pair counts in the radial histogram bins.
Shape: \((N_\mathrm{bins},)\).
- results.rdfnumpy.ndarray
One of
norm="rdf": the radial distribution function \(g_{\alpha\beta}(r)\),norm="density": the single particle density \(n_{\alpha\beta}(r)\), ornorm=None: the raw particle pair count in the radial histogram bins.
Shape: \((N_\mathrm{bins},)\).
- results.coordination_numbersnumpy.ndarray
Coordination numbers \(n_k\). Only available after running
calculate_coordination_numbers().- results.pmfnumpy.ndarray
Potential of mean force \(w_{\alpha\beta}(r)\). Only available after running
calculate_pmf().Shape: \((N_\mathrm{bins},)\).
Reference unit: \(\mathrm{kJ/mol}\).
- results.wavenumbersnumpy.ndarray
Wavenumbers \(q\). Only available after running
calculate_structure_factor().Reference unit: \(\mathrm{Å}^{-1}\).
- results.ssfnumpy.ndarray
Static or partial structure factor. Only available after running
calculate_structure_factor().Shape: Same as results.wavenumbers.
Methods
Calculates the coordination numbers \(n_k\).
Calculates the potential of mean force \(w_{\alpha\beta}(r)\).
Computes the static or partial structure factor \(S_{\alpha\beta}(q)\) using the radial histogram bins \(r\) and the radial distribution function \(g_{\alpha\beta}(r)\) for an isotropic fluid.
Tuple with backends supported by the core library for a given class.
Performs the calculation.
Saves results to a binary or archive file in NumPy format.
- calculate_coordination_numbers(rho: float, *, n_coord_nums: int = 2, threshold: float = 0.1) None[source]¶
Calculates the coordination numbers \(n_k\).
If the radial distribution function \(g_{\alpha\beta}(r)\) does not contain \(k\) local minima, this method will return numpy.nan for the coordination numbers that could not be calculated.
- Parameters:
- rhofloat
Number density \(\rho_\beta\) of species \(\beta\).
Reference unit: \(\mathrm{nm}^{-3}\).
- n_coord_numsint, keyword-only, default:
2 Number of coordination numbers \(n_\mathrm{coord}\) to calculate.
- thresholdfloat, keyword-only, default:
0.1 Minimum \(g_{\alpha\beta}(r)\) value that must be reached before local minima are found.
- calculate_pmf(temperature: float | Quantity | Quantity) None[source]¶
Calculates the potential of mean force \(w_{\alpha\beta}(r)\).
- Parameters:
- temperaturefloat or openmm.unit.Quantity
System temperature \(T\).
Note
If
reduced=Truewas set in theRDFconstructor, temperature should be equal to the energy scale. When the Lennard-Jones potential is used, it generally means that \(T^*=1\), or temperature=1.Reference unit: \(\mathrm{K}\).
- calculate_structure_factor(rho: float, x_a: float = None, x_b: float = None, q: ndarray[float] = None, *, q_lower: float = None, q_upper: float = None, n_q: int = 1000, formalism: str = 'FZ') None[source]¶
Computes the static or partial structure factor \(S_{\alpha\beta}(q)\) using the radial histogram bins \(r\) and the radial distribution function \(g_{\alpha\beta}(r)\) for an isotropic fluid.
- Parameters:
- rhofloat
Bulk number density \(\rho\).
Reference unit: \(\mathrm{Å}^{-3}\).
- x_afloat, default:
1 Number concentration of species \(\alpha\). Required if the two atom groups are not identical.
- x_bfloat, optional
Number concentration of species \(\beta\). Required if the two atom groups are not identical.
- qnumpy.ndarray, optional
Wavenumbers \(q\).
Reference unit: \(\mathrm{Å}^{-1}\).
- q_lowerfloat, keyword-only, optional
Lower bound for the wavenumbers \(q\). Has no effect if q is specified.
Reference unit: \(\mathrm{Å}^{-1}\).
- q_upperfloat, keyword-only, optional
Upper bound for the wavenumbers \(q\). Has no effect if q is specified.
Reference unit: \(\mathrm{Å}^{-1}\).
- n_qint, keyword-only, default:
1_000 Number of wavenumbers \(q\) to generate. Has no effect if q is specified.
- formalism :`str`, keyword-only, default: :code:`”FZ”`
Formalism to use for the partial structure factor. Has no effect if the two atom groups are the same.
Valid values:
"general": A general formalism similar to that of the static structure factor, except the second term is multiplied by \(x_ix_j\)."FZ": Faber–Ziman formalism."AL": Ashcroft–Langreth formalism.
See also
For more information, see
calculate_structure_factor().
- classmethod get_supported_backends()¶
Tuple with backends supported by the core library for a given class. User can pass either one of these values as
backend=...torun()method, or a custom object that hasapplymethod (see documentation forrun()):‘serial’: no parallelization
‘multiprocessing’: parallelization using multiprocessing.Pool
‘dask’: parallelization using dask.delayed.compute(). Requires installation of mdanalysis[dask]
If you want to add your own backend to an existing class, pass a
backends.BackendBasesubclass (see its documentation to learn how to implement it properly), and specifyunsupported_backend=True.- Returns:
- tuple
names of built-in backends that can be used in
run(backend=...)()
Added in version 2.8.0: ..
- property parallelizable¶
Boolean mark showing that a given class can be parallelizable with split-apply-combine procedure. Namely, if we can safely distribute
_single_frame()to multiple workers and then combine them with a proper_conclude()call. If set toFalse, no backends except forserialare supported.Note
If you want to check parallelizability of the whole class, without explicitly creating an instance of the class, see
_analysis_algorithm_is_parallelizable. Note that you setting it to other value will break things if the algorithm behind the analysis is not trivially parallelizable.- Returns:
- bool
if a given
AnalysisBasesubclass instance is parallelizable with split-apply-combine, or not
Added in version 2.8.0: ..
- run(start: int = None, stop: int = None, step: int = None, frames: slice | ndarray[int] = None, verbose: bool = None, **kwargs) SerialAnalysisBase | ParallelAnalysisBase¶
Performs the calculation.
See also
For parallel-specific keyword arguments, see
ParallelAnalysisBase.run().- Parameters:
- startint, optional
Starting frame for analysis.
- stopint, optional
Ending frame for analysis.
- stepint, optional
Number of frames to skip between each analyzed frame.
- framesslice or array-like, optional
Index or logical array of the desired trajectory frames.
- verbosebool, optional
Determines whether detailed progress is shown.
- **kwargs
Additional keyword arguments to pass to
MDAnalysis.lib.log.ProgressBar.
- Returns:
- selfSerialAnalysisBase or ParallelAnalysisBase
Analysis object with results.
- save(file: str | TextIO, archive: bool = True, compress: bool = True, **kwargs) None¶
Saves results to a binary or archive file in NumPy format.
- Parameters:
- filestr or file
Filename or file-like object where the data will be saved. If file is a str, the
.npyor.npzextension will be appended automatically if not already present.- archivebool, default:
True Determines whether the results are saved to a single archive file. If True, the data is stored in a
.npzfile. Otherwise, the data is saved to multiple.npyfiles.- compressbool, default:
True Determines whether the
.npzfile is compressed. Has no effect whenarchive=False.- **kwargs
Additional keyword arguments to pass to
numpy.save(),numpy.savez(), ornumpy.savez_compressed(), depending on the values of archive and compress.