pymatgen.analysis.diffusion.analyzer module¶
A module to perform diffusion analyses.
For example, calculating diffusivity from mean square displacements etc.). If you use this module, please consider citing the following papers:
Ong, S. P., Mo, Y., Richards, W. D., Miara, L., Lee, H. S., & Ceder, G.
(2013). Phase stability, electrochemical stability and ionic conductivity
of the Li10+-1MP2X12 (M = Ge, Si, Sn, Al or P, and X = O, S or Se) family
of superionic conductors. Energy & Environmental Science, 6(1), 148.
doi:10.1039/c2ee23355j
Mo, Y., Ong, S. P., & Ceder, G. (2012). First Principles Study of the
Li10GeP2S12 Lithium Super Ionic Conductor Material. Chemistry of Materials,
24(1), 15-17. doi:10.1021/cm203303y
- class DiffusionAnalyzer(structure, displacements, specie, temperature, time_step, step_skip, smoothed='max', min_obs=30, avg_nsteps=1000, lattices=None, c_ranges=None, c_range_include_edge=False, structures=None)[source]¶
Bases:
MSONable
Class for performing diffusion analysis.
This constructor is meant to be used with pre-processed data. Other convenient constructors are provided as class methods (see from_vaspruns and from_files).
Given a matrix of displacements (see arguments below for expected format), the diffusivity is given by:
D = 1 / 2dt * <mean square displacement>
where d is the dimensionality, t is the time. To obtain a reliable diffusion estimate, a least squares regression of the MSD against time to obtain the slope, which is then related to the diffusivity.
For traditional analysis, use smoothed=False and weighted=False.
- Parameters:
structure (Structure) – Initial structure.
displacements (array) – Numpy array of with shape [site, time step, axis]
specie (Element/Species) – Species to calculate diffusivity for as a String. E.g., “Li”.
temperature (float) – Temperature of the diffusion run in Kelvin.
time_step (int) – Time step between measurements.
step_skip (int) – Sampling frequency of the displacements ( time_step is multiplied by this number to get the real time between measurements)
smoothed (str) –
Whether to smooth the MSD, and what mode to smooth. Supported modes are:
”max”, which tries to use the maximum # of data points for each time origin, subject to a minimum # of observations given by min_obs, and then weights the observations based on the variance accordingly. This is the default.
”constant”, in which each timestep is averaged over the number of time_steps given by min_steps.
None / False / any other false-like quantity. No
smoothing.
min_obs (int) – Used with smoothed=”max”. Minimum number of observations to have before including in the MSD vs dt calculation. E.g. If a structure has 10 diffusing atoms, and min_obs = 30, the MSD vs dt will be calculated up to dt = total_run_time / 3, so that each diffusing atom is measured at least 3 uncorrelated times. Only applies in smoothed=”max”.
avg_nsteps (int) – Used with smoothed=”constant”. Determines the number of time steps to average over to get the msd for each timestep. Default of 1000 is usually pretty good.
lattices (array) – Numpy array of lattice matrix of every step. Used for NPT-AIMD. For NVT-AIMD, the lattice at each time step is set to the lattice in the “structure” argument.
c_ranges (list) – A list of fractional ranges of z-axis to define regions and calculate regional MSD and regional diffusivity. If the start and end positions of a diffusing specie between two time steps are all within the c_ranges, that displacement is collected to calculate MSD in that time step. units: Å. Default to None.
c_range_include_edge (bool) – Whether to include displacements start or end on the edge of the defined c_ranges into the calculation of regional MSD. Default to False.
structures (list) – A list of trajectory structures only used in the calculation of regional MSD and regional diffusivity. These structures should be the same as those used to construct the diffusion analyzer. Default to None.
- export_msdt(filename)[source]¶
Writes MSD data to a csv file that can be easily plotted in other software.
- 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_dict(d)[source]¶
- Parameters:
d (dict) – Dict representation.
Returns: DiffusionAnalyzer
- classmethod from_files(filepaths, specie, step_skip=10, ncores=None, initial_disp=None, initial_structure=None, **kwargs)[source]¶
Convenient constructor that takes in a list of vasprun.xml paths to perform diffusion analysis.
- Parameters:
filepaths ([str]) – List of paths to vasprun.xml files of runs. ( must be ordered in sequence of MD simulation). For example, you may have done sequential VASP runs and they are in run1, run2, run3, etc. You should then pass in [“run1/vasprun.xml”, “run2/vasprun.xml”, …].
specie (Element/Species) – Species to calculate diffusivity for as a String. E.g., “Li”.
step_skip (int) – Sampling frequency of the displacements ( time_step is multiplied by this number to get the real time between measurements)
ncores (int) – Numbers of cores to use for multiprocessing. Can speed up vasprun parsing considerably. Defaults to None, which means serial. It should be noted that if you want to use multiprocessing, the number of ionic steps in all vasprun .xml files should be a multiple of the ionic_step_skip. Otherwise, inconsistent results may arise. Serial mode has no such restrictions.
initial_disp (np.ndarray) – Sometimes, you need to iteratively compute estimates of the diffusivity. This supplies an initial displacement that will be added on to the initial displacements. Note that this makes sense only when smoothed=False.
initial_structure (Structure) – Like initial_disp, this is used for iterative computations of estimates of the diffusivity. You typically need to supply both variables. This stipulates the initial structure from which the current set of displacements are computed.
**kwargs – kwargs supported by the :class:`DiffusionAnalyzer`_. Examples include smoothed, min_obs, avg_nsteps.
- classmethod from_structures(structures, specie, temperature, time_step, step_skip, initial_disp=None, initial_structure=None, **kwargs)[source]¶
Convenient constructor that takes in a list of Structure objects to perform diffusion analysis.
- Parameters:
structures ([Structure]) – list of Structure objects (must be ordered in sequence of run). E.g., you may have performed sequential VASP runs to obtain sufficient statistics.
specie (Element/Species) – Species to calculate diffusivity for as a String. E.g., “Li”.
temperature (float) – Temperature of the diffusion run in Kelvin.
time_step (int) – Time step between measurements.
step_skip (int) – Sampling frequency of the displacements ( time_step is multiplied by this number to get the real time between measurements)
initial_disp (np.ndarray) – Sometimes, you need to iteratively compute estimates of the diffusivity. This supplies an initial displacement that will be added on to the initial displacements. Note that this makes sense only when smoothed=False.
initial_structure (Structure) – Like initial_disp, this is used for iterative computations of estimates of the diffusivity. You typically need to supply both variables. This stipulates the initial structure from which the current set of displacements are computed.
**kwargs – kwargs supported by the :class:`DiffusionAnalyzer`_. Examples include smoothed, min_obs, avg_nsteps.
- classmethod from_vaspruns(vaspruns, specie, initial_disp=None, initial_structure=None, **kwargs)[source]¶
Convenient constructor that takes in a list of Vasprun objects to perform diffusion analysis.
- Parameters:
vaspruns ([Vasprun]) – List of Vaspruns (must be ordered in sequence of MD simulation). E.g., you may have performed sequential VASP runs to obtain sufficient statistics.
specie (Element/Species) – Species to calculate diffusivity for as a String. E.g., “Li”.
initial_disp (np.ndarray) – Sometimes, you need to iteratively compute estimates of the diffusivity. This supplies an initial displacement that will be added on to the initial displacements. Note that this makes sense only when smoothed=False.
initial_structure (Structure) – Like initial_disp, this is used for iterative computations of estimates of the diffusivity. You typically need to supply both variables. This stipulates the initial stricture from which the current set of displacements are computed.
**kwargs – kwargs supported by the :class:`DiffusionAnalyzer`_. Examples include smoothed, min_obs, avg_nsteps.
- get_drift_corrected_structures(start=None, stop=None, step=None)[source]¶
Returns an iterator for the drift-corrected structures. Use of iterator is to reduce memory usage as # of structures in MD can be huge. You don’t often need all the structures all at once.
- Parameters:
start (int) – Applies a start to the iterator.
stop (int) – Applies a stop to the iterator.
step (int) – Applies a step to the iterator. Faster than applying it after generation, as it reduces the number of structures created.
- get_framework_rms_plot(granularity=200, matching_s=None)[source]¶
Get the plot of rms framework displacement vs time. Useful for checking for melting, especially if framework atoms can move via paddle-wheel or similar mechanism (which would show up in max framework displacement but doesn’t constitute melting).
- Parameters:
granularity (int) – Number of structures to match
matching_s (Structure) – Optionally match to a disordered structure instead of the first structure in the analyzer. Required when a secondary mobile ion is present.
Notes
The method doesn’t apply to NPT-AIMD simulation analysis.
- get_msd_plot(mode='specie')[source]¶
Get the plot of the smoothed msd vs time graph. Useful for checking convergence. This can be written to an image file.
- Parameters:
mode (str) – Determines type of msd plot. By “species”, “sites”, or direction (default). If mode = “mscd”, the smoothed mscd vs. time will be plotted.
- get_summary_dict(include_msd_t=False, include_mscd_t=False)[source]¶
Provides a summary of diffusion information.
- Parameters:
include_msd_t (bool) – Whether to include mean square displace and time data with the data.
include_mscd_t (bool) – Whether to include mean square charge displace and time data with the data.
- Returns:
(dict) of diffusion and conductivity data.
- plot_msd(mode='default')[source]¶
Plot the smoothed msd vs time graph. Useful for checking convergence.
- Parameters:
mode (str) – Can be “default” (the default, shows only the MSD for the diffusing specie, and its components), “ions” (individual square displacements of all ions), “species” (mean square displacement by specie), or “mscd” (overall mean square charge displacement for diffusing specie).
- fit_arrhenius(temps: Sequence[float], diffusivities: Sequence[float], mode: Literal['linear', 'exp'] = 'linear', diffusivity_errors: Sequence[float] | None = None) tuple[float, float, float | None] [source]¶
- Returns Ea, c, standard error of Ea from the Arrhenius fit.
D = c * exp(-Ea/kT).
- Parameters:
temps ([float]) – A sequence of temperatures. units: K
diffusivities ([float]) – A sequence of diffusivities (e.g., from DiffusionAnalyzer.diffusivity). units: cm^2/s
mode (str) –
- The fitting mode. Supported modes are:
”linear” (default), which fits ln(D) vs 1/T.
”exp”, which fits D vs T.
Hint: Use “exp” with diffusivity errors if the errors are not homoscadastic in the linear representation. Avoid using “exp” without errors.
diffusivity_errors ([float]) – A sequence of absolute errors in diffusivities. units: cm^2/s
- get_arrhenius_plot(temps, diffusivities, diffusivity_errors=None, mode: Literal['linear', 'exp'] = 'linear', unit: Literal['eV', 'meV'] = 'meV', **kwargs)[source]¶
Returns an Arrhenius plot.
- Parameters:
temps ([float]) – A sequence of temperatures.
diffusivities ([float]) – A sequence of diffusivities (e.g., from DiffusionAnalyzer.diffusivity).
diffusivity_errors ([float]) – A sequence of errors for the diffusivities. If None, no error bar is plotted.
mode (str) – The fitting mode. See fit_arrhenius for details.
unit (str) – The unit for the activation energy. Supported units are “eV” and “meV”.
**kwargs – Any keyword args supported by matplotlib.pyplot.plot.
- Returns:
A matplotlib.Axes object. Do ax.show() to show the plot.
- get_conversion_factor(structure, species, temperature)[source]¶
Conversion factor to convert between cm^2/s diffusivity measurements and mS/cm conductivity measurements based on number of atoms of diffusing species. Note that the charge is based on the oxidation state of the species (where available), or else the number of valence electrons (usually a good guess, esp for main group ions).
- Parameters:
structure (Structure) – Input structure.
species (Element/Species) – Diffusing species.
temperature (float) – Temperature of the diffusion run in Kelvin.
- Returns:
Conversion factor. Conductivity (in mS/cm) = Conversion Factor * Diffusivity (in cm^2/s)
- get_diffusivity_from_msd(msd, dt, smoothed='max')[source]¶
Returns diffusivity and standard deviation of diffusivity.
D = 1 / 2dt * <mean square displacement>
where d is the dimensionality, t is the time. To obtain a reliable diffusion estimate, a least squares regression of the MSD against time to obtain the slope, which is then related to the diffusivity.
For traditional analysis, use smoothed=False.
- Parameters:
msd ([float]) – A sequence of mean square displacements. units: Å^2
dt ([float]) – A sequence of time steps corresponding to MSD. units: fs
smoothed (str) –
Whether to smooth the MSD, and what mode to smooth. Supported modes are:
”max”, which tries to use the maximum # of data points for each time origin, subject to a minimum # of observations given by min_obs, and then weights the observations based on the variance accordingly. This is the default.
”constant”, in which each timestep is averaged over the number of time_steps given by min_steps.
None / False / any other false-like quantity. No
smoothing.
- get_extrapolated_conductivity(temps, diffusivities, new_temp, structure, species)[source]¶
Returns extrapolated mS/cm conductivity.
- Parameters:
temps ([float]) – A sequence of temperatures. units: K
diffusivities ([float]) – A sequence of diffusivities (e.g., from DiffusionAnalyzer.diffusivity). units: cm^2/s
new_temp (float) – desired temperature. units: K
structure (structure) – Structure used for the diffusivity calculation
species (string/Species) – conducting species
- Returns:
(float) Conductivity at extrapolated temp in mS/cm.
- get_extrapolated_diffusivity(temps, diffusivities, new_temp, mode: Literal['linear', 'exp'] = 'linear')[source]¶
Returns (Arrhenius) extrapolated diffusivity at new_temp.
- Parameters:
temps ([float]) – A sequence of temperatures. units: K
diffusivities ([float]) – A sequence of diffusivities (e.g., from DiffusionAnalyzer.diffusivity). units: cm^2/s
new_temp (float) – desired temperature. units: K
mode (str) – The fitting mode. See fit_arrhenius for details.
- Returns:
(float) Diffusivity at extrapolated temp in cm^2/s.