pymatgen.analysis.diffusion.aimd.rdf module

RDF implementation.

class RadialDistributionFunction(structures: list, indices: list, reference_indices: list, ngrid: int = 101, rmax: float = 10.0, cell_range: int = 1, sigma: float = 0.1)[source]

Bases: object

Calculate the average radial distribution function for a given set of structures.

Parameters:
  • structures ([Structure]) – list of structure objects with the same composition. Allow for ensemble averaging.

  • ngrid (int) – Number of radial grid points.

  • rmax (float) – Maximum of radial grid (the minimum is always zero) in Angstrom.

  • cell_range (int) – Range of translational vector elements associated with supercell. Default is 1, i.e. including the adjacent image cells along all three directions.

  • sigma (float) – Smearing of a Gaussian function.

  • indices ([int]) – A list of atom index of interest.

  • reference_indices ([int]) – set this option along with ‘indices’ parameter to compute radial distribution function.

property coordination_number

returns running coordination number.

Returns:

numpy array

export_rdf(filename: str)[source]

Output RDF data to a csv file.

Parameters:

filename (str) – Filename. Supported formats are csv and dat. If the extension is csv, a csv file is written. Otherwise, a dat format is assumed.

classmethod from_species(structures: list, ngrid: int = 101, rmax: float = 10.0, cell_range: int = 1, sigma: float = 0.1, species: tuple | list = ('Li', 'Na'), reference_species: tuple | list | None = None)[source]

Initialize using species.

Parameters:
  • structures (list of pmg_structure objects) – list of structure objects with the same composition. Allow for ensemble averaging.

  • ngrid (int) – Number of radial grid points.

  • rmax (float) – Maximum of radial grid (the minimum is always zero).

  • cell_range (int) – Range of translational vector elements associated with supercell. Default is 1, i.e. including the adjacent image cells along all three directions.

  • sigma (float) – Smearing of a Gaussian function.

  • species (list[string]) – A list of specie symbols of interest.

  • reference_species (list[string]) –

    set this option along with ‘species’ parameter to compute radial distribution function. eg: species=[“H”], reference_species=[“O”] to compute

    O-H pair distribution in a water MD simulation.

get_rdf_plot(label: str | None = None, xlim: tuple = (0.0, 8.0), ylim: tuple = (-0.005, 3.0), loc_peak: bool = False)[source]

Plot the average RDF function.

Parameters:
  • label (str) – The legend label.

  • xlim (list) – Set the x limits of the current axes.

  • ylim (list) – Set the y limits of the current axes.

  • loc_peak (bool) – Label peaks if True.

class RadialDistributionFunctionFast(structures: Structure | list[Structure], rmin: float = 0.0, rmax: float = 10.0, ngrid: float = 101, sigma: float = 0.0, n_jobs=None)[source]

Bases: object

Fast radial distribution analysis.

This method calculates rdf on np.linspace(rmin, rmax, ngrid) points.

Parameters:
  • structures (list of pymatgen Structures) – structures to compute RDF

  • rmin (float) – minimal radius

  • rmax (float) – maximal radius

  • ngrid (int) – number of grid points, defaults to 101

  • sigma (float) – smooth parameter

  • n_jobs (int) – number of CPUs in processing

get_coordination_number(ref_species, species, is_average=True)[source]

returns running coordination number.

Parameters:
  • ref_species (list of species or just single specie str) – the reference species. The rdfs are calculated with these species at the center

  • species (list of species or just single specie str) – the species that we are interested in. The rdfs are calculated on these species.

  • is_average (bool) – whether to take structural average

Returns:

numpy array

get_one_rdf(ref_species: str | list[str], species: str | list[str], index=0)[source]

Get the RDF for one structure, indicated by the index of the structure in all structures.

Parameters:
  • ref_species (list of species or just single specie str) – the reference species. The rdfs are calculated with these species at the center

  • species (list of species or just single specie str) – the species that we are interested in. The rdfs are calculated on these species.

  • index (int) – structure index in the list

Returns:

(x, rdf) x is the radial points, and rdf is the rdf value.

get_rdf(ref_species: str | list[str], species: str | list[str], is_average=True)[source]
Parameters:
  • ref_species (list of species or just single specie str) – the reference species. The rdfs are calculated with these species at the center

  • species (list of species or just single specie str) – the species that we are interested in. The rdfs are calculated on these species.

  • is_average (bool) – whether to take the average over all structures.

Returns:

(x, rdf) x is the radial points, and rdf is the rdf value.