kinisi.analyze#

This module contains the API classes for kinisi. It is anticipated that this is where the majority of interaction with the package will occur. This module includes:

  • the kinisi.analyze.DiffusionAnalyzer class for MSD and diffusion analysis;

  • the kinisi.analyze.JumpDiffusionAnalyzer class for MSTD and collective diffusion analysis; and

  • the kinisi.analyze.ConductivityAnalyzer class for MSCD and conductivity analysis.

These are all compatible with VASP Xdatcar output files, pymatgen structures and any MD trajectory that the MDAnalysis package can handle.

The kinisi.analyze.DiffusionAnalyser class enable the evaluation of tracer mean-squared displacment and the self-diffusion coefficient.

class kinisi.diffusion_analyzer.DiffusionAnalyzer(trajectory)#

Bases: Analyzer

The class for the investigation of the self-diffusion.

Parameters:
  • trajectory (Parser) – The parsed trajectory from some input file. This will be of type Parser, but the specifics depend on the parser that is used.

  • da – The mean-squared displacement data, which is a scipp.DataArray object. This is calculated from the trajectory data and is used to determine the diffusion coefficient.

classmethod from_xdatcar(trajectory, specie, time_step, step_skip, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, masses=None, progress=True)#

Constructs the necessary kinisi objects for analysis from a single or a list of pymatgen.io.vasp.outputs.Xdatcar objects.

Parameters:
  • trajectory (Union[pymatgen.io.vasp.outputs.Xdatcar, list[pymatgen.io.vasp.outputs.Xdatcar]]) – The pymatgen.io.vasp.outputs.Xdatcar or list of these that should be parsed.

  • specie (Union[pymatgen.core.periodic_table.Element, pymatgen.core.periodic_table.Specie]) – Specie to calculate diffusivity for as a String, e.g. 'Li'.

  • time_step (Variable) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (Variable) – Sampling freqency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • dtype (str | None) – If trajectory is a pymatgen.io.vasp.outputs.Xdatcar object, this should be None. However, if a list of pymatgen.io.vasp.outputs.Xdatcar objects is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.:

  • dt (Variable) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • specie_indices (Variable) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (Variable) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

DiffusionAnalyzer

Returns:

The DiffusionAnalyzer object with the mean-squared displacement calculated.

classmethod from_universe(trajectory, specie=None, time_step=None, step_skip=None, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, masses=None, progress=True)#

Constructs the necessary kinisi objects for analysis from a MDAnalysis.Universe object.

Parameters:
  • trajectory (MDAnalysis.core.universe.Universe) – The pymatgen.io.vasp.outputs.Xdatcar or list of these that should be parsed.

  • specie (str) – Specie to calculate diffusivity for as a String, e.g. 'Li'.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling freqency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • dtype (str | None) – If trajectory is a pymatgen.io.vasp.outputs.Xdatcar object, this should be None. However, if a list of pymatgen.io.vasp.outputs.Xdatcar objects is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.:

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • specie_indices (Variable) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (Variable) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

DiffusionAnalyzer

Returns:

The DiffusionAnalyzer object with the mean-squared displacement calculated.

classmethod from_ase(trajectory, specie=None, time_step=None, step_skip=None, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, masses=None, progress=True)#

Constructs the necessary kinisi objects for analysis from an ASE trajectory.

Parameters:
  • trajectory (Union[ase.io.trajectory.Trajectory, list[ase.io.trajectory.Trajectory]]) – The ASE trajectory or list of ASE trajectories to parse.

  • specie (str) – Specie to calculate diffusivity for as a String, e.g. 'Li'.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling frequency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • dtype (str | None) – If trajectory is a single ASE trajectory, this should be None. However, if a list of ASE trajectories is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • specie_indices (Variable) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (Variable) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

DiffusionAnalyzer

Returns:

The DiffusionAnalyzer object with the mean-squared displacement calculated.

diffusion(start_dt, cond_max=1e+16, recondition=False, fit_intercept=True, n_samples=1000, n_walkers=32, n_burn=500, n_thin=10, progress=True, random_state=None)#

Calculate the diffusion coefficient using the mean-squared displacement data.

Parameters:
  • start_dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The time at which the diffusion regime begins.

  • cond_max (float) – The maximum condition number of the covariance matrix. Optional, default is 1e16.

  • recondition (bool) – Whether to recondition the covariance matrix. Optional, default is False.

  • fit_intercept (bool) – Whether to fit an intercept. Optional, default is True.

  • n_samples (int) – The number of MCMC samples to take. Optional, default is 1000.

  • n_walkers (int) – The number of walkers to use in the MCMC. Optional, default is 32.

  • n_burn (int) – The number of burn-in samples to discard. Optional, default is 500.

  • n_thin (int) – The thinning factor for the MCMC samples. Optional, default is 10.

  • progress (bool) – Whether to show the progress bar. Optional, default is True.

  • random_state (RandomState) – The random state to use for the MCMC. Optional, default is None.

Return type:

None

property distributions: array#
Returns:

A distribution of samples for the linear relationship that can be used for easy

plotting of credible intervals.

property D: VariableLikeType#
Returns:

The self-diffusion coefficient.

property msd: VariableLikeType#
Returns:

The mean-squared displacement.

property flatchain: DataGroup#
Returns:

The flatchain of the MCMC samples.

The kinisi.analyze.JumpDiffusionAnalyzer class allows the study of jump diffusion and the collective motion of particles.

class kinisi.jump_diffusion_analyzer.JumpDiffusionAnalyzer(trajectory)#

Bases: Analyzer

The kinisi.analyze.JumpDiffusionAnalyzer class performs analysis of collective diffusion relationships in materials. This is achieved through the application of a bootstrapping methodology to obtain the most statistically accurate values for total mean squared displacement uncertainty and covariance. The time-dependence of the MSTD is then modelled in a generalised least squares fashion to obtain the jump diffusion coefficient and offset using Markov chain Monte Carlo maximum likelihood sampling.

Parameters:
  • trajectory (Parser) – The parsed trajectory from some input file. This will be of type Parser, but the specifics depend on the parser that is used.

  • da – The mean-squared total displacement data, which is a scipp.DataArray object. This is calculated from the trajectory data and is used to determine the jump diffusion coefficient.

classmethod from_xdatcar(trajectory, specie, time_step, step_skip, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, masses=None, system_particles=1, progress=True)#

Constructs the necessary kinisi objects for analysis from a single or a list of pymatgen.io.vasp.outputs.Xdatcar objects.

Parameters:
  • trajectory (Union[pymatgen.io.vasp.outputs.Xdatcar, list[pymatgen.io.vasp.outputs.Xdatcar]]) – The pymatgen.io.vasp.outputs.Xdatcar or list of these that should be parsed.

  • specie (Union[pymatgen.core.periodic_table.Element, pymatgen.core.periodic_table.Specie]) – Specie to calculate diffusivity for as a String, e.g. 'Li'.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling freqency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • dtype (str | None) – If trajectory is a pymatgen.io.vasp.outputs.Xdatcar object, this should be None. However, if a list of pymatgen.io.vasp.outputs.Xdatcar objects is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.:

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • specie_indices (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • system_particles (int) – The number of system particles to average over. Note that the constitution of the system particles are defined in index order, i.e., two system particles will involve splitting the particles down the middle into each. Optional, defaults to 1.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

JumpDiffusionAnalyzer

Returns:

The JumpDiffusionAnalyzer object with the mean-squared total displacement calculated.

classmethod from_universe(trajectory, specie=None, time_step=None, step_skip=None, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, masses=None, system_particles=1, progress=True)#

Constructs the necessary kinisi objects for analysis from a MDAnalysis.core.universe.Universe object.

Parameters:
  • trajectory (MDAnalysis.core.universe.Universe) – The pymatgen.io.vasp.outputs.Xdatcar or list of these that should be parsed.

  • specie (str) – Specie to calculate diffusivity for as a String, e.g. 'Li'.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling freqency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • dtype (str | None) – If trajectory is a pymatgen.io.vasp.outputs.Xdatcar object, this should be None. However, if a list of pymatgen.io.vasp.outputs.Xdatcar objects is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.:

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • specie_indices (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • system_particles (int) – The number of system particles to average over. Note that the constitution of the system particles are defined in index order, i.e., two system particles will involve splitting the particles down the middle into each. Optional, defaults to 1.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

JumpDiffusionAnalyzer

Returns:

The JumpDiffusionAnalyzer object with the mean-squared total displacement calculated.

classmethod from_ase(trajectory, specie=None, time_step=None, step_skip=None, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, masses=None, system_particles=1, progress=True)#

Constructs the necessary kinisi objects for analysis from an ASE trajectory.

Parameters:
  • trajectory (Union[ase.io.trajectory.Trajectory, list[ase.io.trajectory.Trajectory]]) – The ASE trajectory or list of ASE trajectories to parse.

  • specie (str) – Specie to calculate diffusivity for as a String, e.g. 'Li'.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling frequency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • dtype (str | None) – If trajectory is a single ASE trajectory, this should be None. However, if a list of ASE trajectories is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • specie_indices (Variable) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (Variable) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • system_particles (int) – The number of system particles to average over. Note that the constitution of the system particles are defined in index order, i.e., two system particles will involve splitting the particles down the middle into each. Optional, defaults to 1.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

JumpDiffusionAnalyzer

Returns:

The JumpDiffusionAnalyzer object with the mean-squared total displacement calculated.

jump_diffusion(start_dt, cond_max=1e+16, fit_intercept=True, n_samples=1000, n_walkers=32, n_burn=500, n_thin=10, progress=True, random_state=None)#

Calculate the diffusion coefficient using the mean-squared displacement data.

Parameters:
  • start_dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The time at which the diffusion regime begins.

  • cond_max (float) – The maximum condition number of the covariance matrix. Optional, default is 1e16.

  • fit_intercept (bool) – Whether to fit an intercept. Optional, default is True.

  • n_samples (int) – The number of MCMC samples to take. Optional, default is 1000.

  • n_walkers (int) – The number of walkers to use in the MCMC. Optional, default is 32.

  • n_burn (int) – The number of burn-in samples to discard. Optional, default is 500.

  • n_thin (int) – The thinning factor for the MCMC samples. Optional, default is 10.

  • progress (bool) – Whether to show the progress bar. Optional, default is True.

  • random_state (RandomState) – The random state to use for the MCMC. Optional, default is None.

Return type:

None

property distributions: array#
Returns:

A distribution of samples for the linear relationship that can be used for easy

plotting of credible intervals.

property D_J: VariableLikeType#
Returns:

The jump diffusion coefficient.

property mstd: VariableLikeType#
Returns:

The mean-squared total displacement.

property flatchain: DataGroup#
Returns:

The flatchain of the MCMC samples.

The :py:class:kinisi.analyze.ConductivityAnalyzer` class for MSCD and conductivity analysis.

class kinisi.conductivity_analyzer.ConductivityAnalyzer(trajectory)#

Bases: Analyzer

The kinisi.analyze.ConductivityAnalyzer class performs analysis of conductivity in materials. This is achieved through the application of a Bayesian regression methodology to obtain optimal estimates of the mean-squared charge displacement and covariance. This is then sampled with a linear model to estimate the conductivity, via the jump diffusion coefficient.

Parameters:

trajectory (Parser) – The parsed trajectory from some input file. This will be of type Parser, but the specifics depend on the parser that is used.

classmethod from_xdatcar(trajectory, specie, time_step, step_skip, ionic_charge, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), species_indices=None, masses=None, system_particles=1, progress=True)#

Constructs the necessary kinisi objects for analysis from a single or a list of pymatgen.io.vasp.outputs.Xdatcar objects.

Parameters:
  • trajectory (Union[pymatgen.io.vasp.outputs.Xdatcar, list[pymatgen.io.vasp.outputs.Xdatcar]]) – The pymatgen.io.vasp.outputs.Xdatcar or list of these that should be parsed.

  • specie (Union[pymatgen.core.periodic_table.Element, pymatgen.core.periodic_table.Specie]) – Specie to calculate diffusivity for as a String, e.g. 'Li'.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling freqency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • ionic_charge (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The ionic charge of the species of interest. This should be either a scipp scalar if all of the ions have the same charge or an array of the charge for each indiviudal ion.

  • dtype (str | None) – If trajectory is a pymatgen.io.vasp.outputs.Xdatcar object, this should be None. However, if a list of pymatgen.io.vasp.outputs.Xdatcar objects is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.:

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • specie_indices – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • system_particles (int) – The number of system particles to average over. Note that the constitution of the system particles are defined in index order, i.e., two system particles will involve splitting the particles down the middle into each. Optional, defaults to 1.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

ConductivityAnalyzer

Returns:

The ConductivityAnalyzer object with the mean-squared charge displacement calculated.

classmethod from_ase(trajectory, specie=None, time_step=None, step_skip=None, ionic_charge=None, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), species_indices=None, masses=None, system_particles=1, progress=True)#

Constructs the necessary kinisi objects for analysis from an ASE trajectory.

Parameters:
  • trajectory (Union[ase.io.trajectory.Trajectory, list[ase.io.trajectory.Trajectory]]) – The ASE trajectory or list of ASE trajectories to parse.

  • specie (Union[str, ase.Atom, None]) – Specie to calculate conductivity for as a String, e.g. 'Li', or an ase.Atom object.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling frequency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • ionic_charge (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The ionic charge of the species of interest. This should be either a scipp scalar if all of the ions have the same charge or an array of the charge for each individual ion.

  • dtype (str | None) – If trajectory is a single ASE trajectory, this should be None. However, if a list of ASE trajectories is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • species_indices (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • system_particles (int) – The number of system particles to average over. Note that the constitution of the system particles are defined in index order, i.e., two system particles will involve splitting the particles down the middle into each. Optional, defaults to 1.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

ConductivityAnalyzer

Returns:

The ConductivityAnalyzer object with the mean-squared charge displacement calculated.

classmethod from_universe(trajectory, specie=None, time_step=None, step_skip=None, ionic_charge=None, dtype=None, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), species_indices=None, masses=None, system_particles=1, progress=True)#

Constructs the necessary kinisi objects for analysis from a MDAnalysis.Universe object.

Parameters:
  • trajectory (MDAnalysis.core.universe.Universe) – The MDAnalysis.Universe object that should be parsed.

  • specie (str) – Specie to calculate conductivity for as a String, e.g. 'Li'.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The input simulation time step, i.e., the time step for the molecular dynamics integrator. Note, that this must be given as a scipp-type scalar. The unit used for the time_step, will be the unit that is use for the time interval values.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling frequency of the simulation trajectory, i.e., how many time steps exist between the output of the positions in the trajectory. Similar to the time_step, this parameter must be a scipp scalar. The units for this scalar should be dimensionless.

  • ionic_charge (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The ionic charge of the species of interest. This should be either a scipp scalar if all of the ions have the same charge or an array of the charge for each individual ion.

  • dtype (str | None) – If trajectory is a MDAnalysis.Universe object, this should be None. However, if a list of MDAnalysis.Universe objects is passed, then it is necessary to identify if these constitute a series of consecutive trajectories or a series of identical starting points with different random seeds, in which case the dtype should be either consecutive or identical.

  • dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time intervals to calculate the displacements over. Optional, defaults to a scipp array ranging from the smallest interval (i.e., time_step * step_skip) to the full simulation length, with a step size the same as the smallest interval.

  • dimension (str) – Dimension/s to find the displacement along, this should be some subset of ‘xyz’ indicating the axes of interest. Optional, defaults to ‘xyz’.

  • distance_unit (Unit) – The unit of distance in the simulation input. This should be a scipp unit and defaults to sc.units.angstrom.

  • species_indices (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The indices of the species to compute the centre of mass of. Optional, defaults to None, which means that all species are considered.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The masses of the species to calculate the diffusion for. Optional, defaults to None, which means that the masses are not considered.

  • system_particles (int) – The number of system particles to average over. Note that the constitution of the system particles are defined in index order, i.e., two system particles will involve splitting the particles down the middle into each. Optional, defaults to 1.

  • progress (bool) – Print progress bars to screen. Optional, defaults to True.

Return type:

ConductivityAnalyzer

Returns:

The ConductivityAnalyzer object with the mean-squared charge displacement calculated.

conductivity(start_dt, temperature, cond_max=1e+16, fit_intercept=True, n_samples=1000, n_walkers=32, n_burn=500, n_thin=10, progress=True, random_state=None)#

Calculation of the conductivity. Keyword arguments will be passed of the bayesian_regression() method.

Parameters:
  • start_dt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The time at which the diffusion regime begins.

  • temperature (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The temperature of the system.

  • cond_max (float) – The maximum condition number of the covariance matrix. Optional, default is 1e16.

  • fit_intercept (bool) – Whether to fit an intercept. Optional, default is True.

  • n_samples (int) – The number of MCMC samples to take. Optional, default is 1000.

  • n_walkers (int) – The number of walkers to use in the MCMC. Optional, default is 32.

  • n_burn (int) – The number of burn-in samples to discard. Optional, default is 500.

  • n_thin (int) – The thinning factor for the MCMC samples. Optional, default is 10.

  • progress (bool) – Whether to show the progress bar. Optional, default is True.

  • random_state (RandomState) – The random state to use for the MCMC. Optional, default is None.

Return type:

None

property distributions: array#
Returns:

A distribution of samples for the linear relationship that can be used for easy

plotting of credible intervals.

property sigma: VariableLikeType#
Returns:

The conductivity.

property mscd: VariableLikeType#
Returns:

The mean-squared charge displacement.

property flatchain: DataGroup#
Returns:

The flatchain of the MCMC samples.

The Analyzer class is the main API interface for users to kinisi. Typically, it is expected that the user would access a specific Analyzer, such as the DiffusionAnalyzer, which sub-class this class.

Note, that all of the classmethod functions for the Analyzer are intended for internal use.

class kinisi.analyzer.Analyzer(trajectory)#

Bases: object

This class is the superclass for the kinisi.analyze.DiffusionAnalyzer, kinisi.analyze.JumpDiffusionAnalyzer and kinisi.analyze.ConductivityAnalyzer classes.

Parameters:

trajectory (Parser) – The parsed trajectory from some input file. This will be of type Parser, but the specifics depend on the parser that is used.

to_hdf5(filename)#

Save the Analyzer object to an HDF5 file.

Parameters:

filename (str) – The name of the file to save the object to.

Return type:

None

classmethod from_hdf5(filename)#

Load the Analyzer object from an HDF5 file.

Parameters:

filename (str) – The name of the file to load the object from.

Return type:

Analyzer

posterior_predictive(n_posterior_samples=None, n_predictive_samples=256, progress=True)#

Sample the posterior predictive distribution. The shape of the resulting array will be (n_posterior_samples * n_predictive_samples, start_dt).

Params posterior_predictive_params:

Parameters for the diffusion.posterior_predictive() method. See the appropriate documentation for more guidence on this dictionary.

Returns:

Samples from the posterior predictive distribution.

property n_atoms: int#
Returns:

The number of atoms in the trajectory.

property intercept: Variable#
Returns:

The intercept of the linear relationship.

property dt: Variable#
Returns:

The time intervals used for the mean-squared displacement.

property da: DataArray#
Returns:

The mean-squared displacement data array.

property dg: DataGroup#
Returns:

The data group containing the analysis information.