pymatgen.analysis.diffusion.aimd.van_hove module¶
Van Hove analysis for correlations.
- class EvolutionAnalyzer(structures: list, rmax: float = 10, step: int = 1, time_step: int = 2)[source]¶
Bases:
object
Analyze the evolution of structures during AIMD simulations.
Initialization the EvolutionAnalyzer from MD simulations. From the structures obtained from MD simulations, we can analyze the structure evolution with time by some quantitative characterization such as RDF and atomic distribution.
If you use this class, please consider citing the following paper: Tang, H.; Deng, Z.; Lin, Z.; Wang, Z.; Chu, I-H; Chen, C.; Zhu, Z.; Zheng, C.; Ong, S. P. “Probing Solid-Solid Interfacial Reactions in All-Solid-State Sodium-Ion Batteries with First-Principles Calculations”, Chem. Mater. (2018), 30(1), pp 163-173.
- Parameters:
structures ([Structure]) – The list of structures obtained from MD simulations.
rmax (float) – Maximum of radial grid (the minimum is always zero).
step (int) – indicate the interval of input structures, which is used to calculated correct time step.
time_step (int) – the time step in fs, POTIM in INCAR.
- static atom_dist(structure: Structure, specie: str, ngrid: int = 101, window: float = 1, direction: str = 'c')[source]¶
Get atomic distribution for a given specie.
- Parameters:
structure (Structure) – input structure
specie (str) – species string for an element
ngrid (int) – Number of radial grid points.
window (float) – number of atoms will be counted within the range (i-window, i+window), unit is angstrom.
direction (str) – Choose from “a”, “b” and “c”. Default is “c”.
- Returns:
atomic concentration along one direction.
- Return type:
density (np.array)
- get_df(func: Callable, save_csv: str | None = None, **kwargs)[source]¶
Get the data frame for a given pair. This step would be very slow if there are hundreds or more structures to parse.
- Parameters:
func (FunctionType) –
structure to spectrum function. choose from rdf (to get radial distribution function, pair required) or get_atomic_distribution (to get atomic distribution, specie required). Extra parameters can be parsed using kwargs. e.g. To get rdf dataframe:
- df = EvolutionAnalyzer.get_df(
func=EvolutionAnalyzer.rdf, pair=(“Na”, “Na”))
- e.g. To get atomic distribution:
- df = EvolutionAnalyzer.get_df(
func=EvolutionAnalyzer.atom_dist, specie=”Na”)
save_csv (str) – save pandas DataFrame to csv.
**kwargs – Pass-through to func.
- Returns:
- index is the radial distance in Angstrom,
and column is the time step in ps.
- Return type:
pandas.DataFrame object
- static get_min_dist(df: DataFrame, tol: float = 1e-10)[source]¶
Get the shortest pair distance from the given DataFrame.
- Parameters:
df (DataFrame) – index is the radial distance in Angstrom, and column is the time step in ps.
tol (float) – any float number less than tol is considered as zero.
- Returns:
The shortest pair distance throughout the table.
- static get_pairs(structure: Structure)[source]¶
Get all element pairs in a structure.
- Parameters:
structure (Structure) – structure
- Returns:
list of tuples
- plot_atomic_evolution(specie: str, direction: str = 'c', cmap=<matplotlib.colors.LinearSegmentedColormap object>, df: ~pandas.core.frame.DataFrame = None)[source]¶
Plot the atomic distribution evolution with time for a given species.
- Parameters:
specie (str) – species string for an element.
direction (str) – Choose from “a”, “b”, “c”. Default to “c”.
cmap (color map) – the color map used in heat map.
df (DataFrame) – external data, index is the atomic distance in Angstrom, and column is the time step in ps.
- Returns:
matplotlib.axes._subplots.AxesSubplot object
- static plot_evolution_from_data(df: ~pandas.core.frame.DataFrame, x_label: str | None = None, cb_label: str | None = None, cmap=<matplotlib.colors.ListedColormap object>)[source]¶
Plot the evolution with time for a given DataFrame. It can be RDF, atomic distribution or other characterization data we might implement in the future.
- Parameters:
df (pandas.DataFrame) – input DataFrame object, index is the radial distance in Angstrom, and column is the time step in ps.
x_label (str) – x label
cb_label (str) – color bar label
cmap (color map) – the color map used in heat map. cmocean.cm.thermal is recommended
- Returns:
matplotlib.axes._subplots.AxesSubplot object
- plot_rdf_evolution(pair: tuple, cmap=<matplotlib.colors.ListedColormap object>, df: ~pandas.core.frame.DataFrame = None)[source]¶
Plot the RDF evolution with time for a given pair.
- Parameters:
pair (str tuple) – e.g. (“Na”, “Na”)
cmap (color map) – the color map used in heat map. cmocean.cm.thermal is recommended
df (DataFrame) – external data, index is the radial distance in Angstrom, and column is the time step in ps.
- Returns:
matplotlib.axes._subplots.AxesSubplot object
- static rdf(structure: Structure, pair: tuple, ngrid: int = 101, rmax: float = 10)[source]¶
Process rdf from a given structure and pair.
- Parameters:
structure (Structure) – input structure.
pair (str tuple) – e.g. (“Na”, “Na”).
ngrid (int) – Number of radial grid points.
rmax (float) – Maximum of radial grid (the minimum is always zero).
- Returns:
rdf (np.array)
- class VanHoveAnalysis(diffusion_analyzer: DiffusionAnalyzer, avg_nsteps: int = 50, ngrid: int = 101, rmax: float = 10.0, step_skip: int = 50, sigma: float = 0.1, cell_range: int = 1, species: tuple | list = ('Li', 'Na'), reference_species: tuple | list | None = None, indices: list | None = None)[source]¶
Bases:
object
Class for van Hove function analysis.
In particular, self-part (Gs) and distinct-part (Gd) of the van Hove correlation function G(r,t) for given species and given structure are computed. If you use this class, please consider citing the following paper:
Zhu, Z.; Chu, I.-H.; Deng, Z. and Ong, S. P. "Role of Na+ Interstitials and Dopants in Enhancing the Na+ Conductivity of the Cubic Na3PS4 Superionic Conductor". Chem. Mater. (2015), 27, pp 8318-8325
Initiation.
- Parameters:
diffusion_analyzer (DiffusionAnalyzer) – A pymatgen.analysis.diffusion.analyzer.DiffusionAnalyzer object
avg_nsteps (int) – Number of t0 used for statistical average
ngrid (int) – Number of radial grid points
rmax (float) – Maximum of radial grid (the minimum is always set zero)
step_skip (int) – # of time steps skipped during analysis. It defines the resolution of the reduced time grid
sigma (float) – Smearing of a Gaussian function
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.
species ([string]) – a list of specie symbols of interest.
reference_species ([string]) – Set this option along with ‘species’ parameter to calculate the distinct-part of van Hove function. Note that the self-part of van Hove function is always computed only for those in “species” parameter.
indices (list of int) – If not None, only a subset of atomic indices will be selected for the analysis. If this is given, “species” parameter will be ignored.
- get_1d_plot(mode: str = 'distinct', times: list | None = None, colors: list | None = None)[source]¶
Plot the van Hove function at given r or t.
- Parameters:
mode (str) – Specify which part of van Hove function to be plotted.
times (list of float) – Time moments (in ps) in which the van Hove function will be plotted.
colors (list strings/tuples) – Additional color settings. If not set, seaborn.color_plaette(“Set1”, 10) will be used.
- Returns:
axes object.
- Return type:
matplotlib.axes._subplots.Axes
- get_3d_plot(figsize: tuple = (12, 8), mode: str = 'distinct')[source]¶
Plot 3D self-part or distinct-part of van Hove function, which is specified by the input argument ‘type’.
- Parameters:
figsize (tuple) – fig size in inches.
mode (str) – ‘distinct’ or ‘both’.
- Returns:
axes object.
- Return type:
matplotlib.axes._subplots.Axes