Diffusion coefficient of molecules using center of mass#
Kinisi
includes the ability to calculate the mean-squared displacement and diffusion coefficient of the center of mass (or geometry) of molecules. This can be done for a single molecule or a collection of molecules. It is important to note that inclusion of rotational motion in the calcuation of diffusion coeffiencents can lead to erronious results. This rotation can be elminated from the calculation by taking the center of mass for each molecule.
[1]:
import numpy as np
from ase.io import read
import matplotlib.pyplot as plt
from kinisi.analyze import DiffusionAnalyzer
We will use a simulation of ethene in ZSM-5 zeolite. This was run in DL_POLY, so we will use ASE
to load in the trajectory (HISTORY) file.
[2]:
traj = read('ethene_zeo_HISTORY.gz', format='dlp-history', index=':')
We want to calculate the diffusion of the center of mass of the ethene molecule. This can be done by setting specie
to None and specifying the indices of the molecules of interest in specie_indices
. To define molecules, a list of lists should be passed under the specie_indices
keyword. The outer list has one entry per molecule and each inter list has the indices of that molecule. Only identical molecules are supported. The masses of the atoms in the molecules can be specified with
masses
. This must be a list with the same length as a molecule (the length of one of the inner lists in specie_indices
).
[3]:
molecules = [[289, 290, 291, 292, 293, 294],
[285, 296, 297, 298, 299, 300]]
masses = [12, 12, 1.008, 1.008, 1.008, 1.008]
p_parms = {'specie': None,
'time_step': 1.2e-03,
'step_skip': 100,
'specie_indices': molecules,
'masses': masses,
'progress': False
}
u_params = {'progress': False}
With the parameters set, we now calcuate the mean squared-displacement.
[4]:
diff = DiffusionAnalyzer.from_ase(traj, parser_params=p_parms, uncertainty_params=u_params)
[5]:
plt.errorbar(diff.dt, diff.msd, diff.msd_std)
plt.ylabel('MSD/Å$^2$')
plt.xlabel(r'$\Delta t$/ps')
plt.show()
[6]:
diff.diffusion(50, {'progress': False})
[7]:
diff.D.n, diff.D.ci()
[7]:
(np.float64(0.00023699364727608642), array([2.09532152e-05, 5.13703213e-04]))
[8]:
diff.intercept.n, diff.intercept.ci()
[8]:
(np.float64(-36.78023639824387), array([-774.80996767, 594.37474368]))
[9]:
diff.D.samples
[9]:
array([2.37446025e-04, 2.87667422e-04, 8.08208437e-05, ...,
4.04911849e-04, 3.19524841e-04, 1.62189867e-04], shape=(3200,))
[10]:
credible_intervals = [[16, 84], [2.5, 97.5], [0.15, 99.85]]
alpha = [0.6, 0.4, 0.2]
plt.plot(diff.dt, diff.msd, 'k-')
for i, ci in enumerate(credible_intervals):
plt.fill_between(diff.dt,
*np.percentile(diff.distribution, ci, axis=1),
alpha=alpha[i],
color='#0173B2',
lw=0)
plt.ylabel('MSD/Å$^2$')
plt.xlabel(r'$\Delta t$/ps')
plt.show()
[11]:
plt.hist(diff.D.samples, density=True)
plt.axvline(diff.D.n, c='k')
plt.xlabel('$D$/cm$^2$s$^{-1}$')
plt.ylabel('$p(D$/cm$^2$s$^{-1})$')
plt.show()