4.8.2.1. Radial Distribution Functions — MDAnalysis.analysis.rdf
This module contains two classes to calculate radial pair distribution functions (radial distribution functions or “RDF”). The RDF \(g_{ab}(r)\) between types of particles \(a\) and \(b\) is
which is normalized so that the RDF becomes 1 for large separations in a homogenous system. The RDF effectively counts the average number of \(b\) neighbours in a shell at distance \(r\) around a \(a\) particle and represents it as a density.
The radial cumulative distribution function is
and the average number of \(b\) particles within radius \(r\)
(with the appropriate density \(\rho\)). The latter function can be used to compute, for instance, coordination numbers such as the number of neighbors in the first solvation shell \(N(r_1)\) where \(r_1\) is the position of the first minimum in \(g(r)\).
We provide options for calculating the density of particle \(b\) in a shell at distance \(r\) around a \(a\) particle, which is
- class MDAnalysis.analysis.rdf.InterRDF(g1, g2, nbins=75, range=(0.0, 15.0), norm='rdf', exclusion_block=None, exclude_same=None, **kwargs)[source]
- Radial distribution function - InterRDFis a tool to calculate average radial distribution functions between two groups of atoms. Suppose we have two AtomGroups- Aand- B.- Acontains atom- A1,- A2, and- Bcontains- B1,- B2. Given- Aand- Bto- InterRDF, the output will be the average of RDFs between- A1and- B1,- A1and- B2,- A2and- B1,- A2and- B2. A typical application is to calculate the RDF of solvent with itself or with another solute.- The radial distribution function is calculated by histogramming distances between all particles in g1 and g2 while taking periodic boundary conditions into account via the minimum image convention. - The exclusion_block keyword may be used to exclude a set of distances from the calculations. - Results are available in the attributes - results.rdfand- results.count.- Parameters:
- g1 (AtomGroup) – First AtomGroup 
- g2 (AtomGroup) – Second AtomGroup 
- nbins (int) – Number of bins in the histogram 
- norm (str, {'rdf', 'density', 'none'}) – - For ‘rdf’ calculate \(g_{ab}(r)\). For ‘density’ the single particle density \(n_{ab}(r)\) is computed. ‘none’ computes the number of particles occurences in each spherical shell. - Added in version 2.3.0. 
- exclusion_block (tuple) – A tuple representing the tile to exclude from the distance array. 
- exclude_same (str) – Will exclude pairs of atoms that share the same “residue”, “segment”, or “chain”. Those are the only valid values. This is intended to remove atoms that are spatially correlated due to direct bonded connections. 
- verbose (bool) – Show detailed progress of the calculation if set to True 
 
 - results.bins
- numpy.ndarrayof the centers of the nbins histogram bins.- Added in version 2.0.0. - Type:
 
 - bins
- Alias to the - results.binsattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.binsinstead.- Type:
 
 - results.edges
- numpy.ndarrayof the nbins + 1 edges of the histogram bins.- Added in version 2.0.0. - Type:
 
 - edges
- Alias to the - results.edgesattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.edgesinstead.- Type:
 
 - results.rdf
- numpy.ndarrayof the radial distribution function values for the- results.bins.- Added in version 2.0.0. - Type:
 
 - rdf
- Alias to the - results.rdfattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.rdfinstead.- Type:
 
 - results.count
- numpy.ndarrayrepresenting the radial histogram, i.e., the raw counts, for all- results.bins.- Added in version 2.0.0. - Type:
 
 - count
- Alias to the - results.countattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.countinstead.- Type:
 
 - Example - First create the - InterRDFobject, by supplying two AtomGroups then use the- run()method- rdf = InterRDF(ag1, ag2) rdf.run() - Results are available through the - results.binsand- results.rdfattributes:- plt.plot(rdf.results.bins, rdf.results.rdf) - The exclusion_block keyword allows the masking of pairs from within the same molecule. For example, if there are 7 of each atom in each molecule, the exclusion mask - (7, 7)can be used.- Added in version 0.13.0. - Changed in version 1.0.0: Support for the start, stop, and step keywords has been removed. These should instead be passed to - InterRDF.run().- Changed in version 2.0.0: Store results as attributes bins, edges, rdf and count of the results attribute of - AnalysisBase.
- class MDAnalysis.analysis.rdf.InterRDF_s(u, ags, nbins=75, range=(0.0, 15.0), norm='rdf', density=False, **kwargs)[source]
- Site-specific radial distribution function - Calculates site-specific radial distribution functions. Instead of two groups of atoms it takes as input a list of pairs of AtomGroup, - [[A, B], [C, D], ...]. Given the same- Aand- Bto- InterRDF_s, the output will be a list of individual RDFs between- A1and- B1,- A1and- B2,- A2and- B1,- A2and- B2(and similarly for- Cand- D). These site-specific radial distribution functions are typically calculated if one is interested in the solvation shells of a ligand in a binding site or the solvation of specific residues in a protein.- Parameters:
- u (Universe) – - a Universe that contains atoms in ags - Deprecated since version 2.3.0: This parameter is superflous and will be removed in MDAnalysis 3.0.0. 
- nbins (int) – Number of bins in the histogram 
- norm (str, {'rdf', 'density', 'none'}) – - For ‘rdf’ calculate \(g_{ab}(r)\). For ‘density’ the single particle density \(n_{ab}(r)\) is computed. ‘none’ computes the number of particles occurences in each spherical shell. - Added in version 2.3.0. 
- density (bool) – - False: calculate \(g_{ab}(r)\); True: calculate the true single particle density \(n_{ab}(r)\). density overwrites the norm parameter. - Added in version 1.0.1: This keyword was available since 0.19.0 but was not documented. Furthermore, it had the opposite meaning. Since 1.0.1 it is officially supported as documented. - Deprecated since version 2.3.0: Instead of density=True use norm=’density’ 
 
 - results.bins
- numpy.ndarrayof the centers of the nbins histogram bins; all individual site-specific RDFs have the same bins.- Added in version 2.0.0. - Type:
 
 - bins
- Alias to the - results.binsattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.binsinstead.- Type:
 
 - results.edges
- array of the - nbins + 1edges of the histogram bins; all individual site-specific RDFs have the same bins.- Added in version 2.0.0. - Type:
 
 - edges
- Alias to the - results.edgesattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.edgesinstead.- Type:
 
 - results.rdf
- listof the site-specific radial distribution functions if norm=’rdf’ or density functions for the- binsif norm=’density’. The list contains- len(ags)entries. Each entry for the- i-th pair [A, B] = ags[i] in ags is a- numpy.ndarraywith shape- (len(A), len(B)), i.e., a stack of RDFs. For example,- results.rdf[i][0, 2]is the RDF between atoms- A[0]and- B[2].- Added in version 2.0.0. - Type:
 
 - rdf
- Alias to the - results.rdfattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.rdfinstead.- Type:
 
 - results.count
- listof the site-specific radial histograms, i.e., the raw counts, for all- results.bins. The data have the same structure as- results.rdfexcept that the arrays contain the raw counts.- Added in version 2.0.0. - Type:
 
 - count
- Alias to the - results.countattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.countinstead.- Type:
 
 - results.cdf
- listof the site-specific cumulative counts, for all- results.bins. The data have the same structure as- results.rdfexcept that the arrays contain the cumulative counts.- This attribute only exists after - get_cdf()has been run.- Added in version 2.0.0. - Type:
 
 - cdf
- Alias to the - results.cdfattribute.- Deprecated since version 2.0.0: This attribute will be removed in 3.0.0. Use - results.cdfinstead.- Type:
 
 - Example - First create the - InterRDF_sobject, by supplying one Universe and one list of pairs of AtomGroups, then use the- run()method:- from MDAnalysisTests.datafiles import GRO_MEMPROT, XTC_MEMPROT u = mda.Universe(GRO_MEMPROT, XTC_MEMPROT) s1 = u.select_atoms('name ZND and resid 289') s2 = u.select_atoms('(name OD1 or name OD2) and resid 51 and sphzone 5.0 (resid 289)') s3 = u.select_atoms('name ZND and (resid 291 or resid 292)') s4 = u.select_atoms('(name OD1 or name OD2) and sphzone 5.0 (resid 291)') ags = [[s1, s2], [s3, s4]] rdf = InterRDF_s(u, ags) rdf.run() - Results are available through the - results.binsand- results.rdfattributes:- plt.plot(rdf.results.bins, rdf.results.rdf[0][0, 0]) - (Which plots the rdf between the first atom in - s1and the first atom in- s2)- To generate the cumulative distribution function (cdf) in the sense of “particles within radius \(r\)”, i.e., \(N_{ab}(r)\), use the - get_cdf()method- cdf = rdf.get_cdf() - Results are available through the - results.cdfattribute:- plt.plot(rdf.results.bins, rdf.results.cdf[0][0, 0]) - (Which plots the cdf between the first atom in - s1and the first atom in- s2)- Added in version 0.19.0. - Changed in version 1.0.0: Support for the start, stop, and step keywords has been removed. These should instead be passed to - InterRDF_s.run().- Changed in version 2.0.0: Store results as attributes bins, edges, rdf, count and cdf of the results attribute of - AnalysisBase.- Changed in version 2.3.0: Introduce norm and exclusion_blocks attributes. - Deprecated since version 2.3.0: Instead of density=True use norm=’density’ - Deprecated since version 2.3.0: The universe parameter is superflous. - get_cdf()[source]
- Calculate the cumulative counts for all sites. - This is the cumulative count within a given radius, i.e., \(N_{ab}(r)\). - The result is returned and also stored in the attribute - results.cdf.- Returns:
- cdf – list of arrays with the same structure as - results.rdf
- Return type: