kinisi.parser#

Parsers for kinisi. This module is responsible for reading in input files from pymatgen, MDAnalysis, and ase.

class kinisi.parser.Parser(coords, latt, time_step, step_skip, dt=None, specie_indices=None, drift_indices=None, masses=None, dimension='xyz')#

Bases: object

The base class for object parsing. :param structure: a pymatgen.core.structure.Structure or a MDAnalysis.core.universe.Universe :type coords: TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]) :param coords: a scipp array with dimensions of time, atom, and dimension), :param lattice: a scipp array with dimensions time,`dimension1`, and dimension2 :param specie: Specie to calculate diffusivity for as a String, e.g. 'Li'. :type time_step: TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]) :param time_step: 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.

Parameters:
  • 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.

  • 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.

  • specie_indices (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Indices of the specie to calculate the diffusivity for. Optional, defaults to None.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Masses of the atoms in the structure. Optional, defaults to None. If used should be a 1D scipp array of dimension ‘atoms in particle’.

  • 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’.

  • progress – Whether to show a progress bar when reading in the structures. Optional, defaults to True.

create_integer_dt(coords, time_step, step_skip)#

Create an integer time interval from the given time intervals (and if necessary the time interval object). Also checks that the time intervals provided in the dt parameter are a valid subset of the simulation time intervals.

Parameters:
  • coords (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The fractional coordiates of the atoms in the trajectory. This should be a scipp array type object with dimensions of ‘particle’, ‘time’, and ‘dimension’.

  • 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.

Raises:

ValueError – If the time intervals provided in the dt parameter are not a subset of the time intervals present in the simulation, based on the time_step and step_skip parameters and number of snapshots in the trajectory.

Return type:

TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])

Returns:

The integer time intervals as a scipp array with dimensions of ‘time interval’.

static orthorhombic_calculate_displacements(coords, lattice)#

Calculate the absolute displacements of the atoms in the trajectory, when the cell is orthorhombic on all frames.

Parameters:
  • coords (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The fractional coordiates of the atoms in the trajectory. This should be a scipp array type object with dimensions of ‘particle’, ‘time’, and ‘dimension’.

  • lattice (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – A series of matrices that describe the lattice in each step in the trajectory. A scipp array with dimensions of ‘time’, ‘dimension1’, and ‘dimension2’.

Return type:

TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])

Returns:

The absolute displacements of the atoms in the trajectory.

static non_orthorhombic_calculate_displacements(coords, lattice)#
Calculate the absolute displacements of the atoms in the trajectory, when a non-orthrhombic cell is used.

This is done by finding the minimum cartesian displacement vector, from its 8 periodic images. This ensures that triclinic cells are treated correctly.

Parameters:
  • coords (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – The fractional coordiates of the atoms in the trajectory. This should be a scipp array type object with dimensions of ‘particle’, ‘time’, and ‘dimension’.

  • lattice (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – A series of matrices that describe the lattice in each step in the trajectory. A scipp array with dimensions of ‘time’, ‘row’, and ‘column’.

Return type:

TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])

Returns:

The absolute displacements of the atoms in the trajectory.

correct_drift(disp)#

Perform drift correction, such that the displacement is calculated normalised to any framework drift.

Parameters:

disp (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Displacements for all atoms in the simulation. A scipp array with dimensions of obs, atom and dimension.

Return type:

TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])

Returns:

Displacements corrected to account for drift of a framework.

abstractmethod get_indices(structure, specie)#
abstractmethod get_drift_indices(structure, specie_indices)#
get_specie_and_drift_indices(specie, specie_indices, drift_indices, structure)#
property coords#

Coordinates of ‘atoms’, this may be the raw coordinates parsed or centres of mass/geometry.

property disp#

Atom displacements, without drift correction.

kinisi.parser.get_molecules(coords, indices, masses)#

Determine framework and non-framework indices for an ase or pymatgen or MDAnalysis compatible file when specie_indices are provided and contain multiple molecules. Warning: This function changes the coords without renaming the object.

Parameters:

of the diffusion. :type masses: TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]) :param masses: Masses associated with indices in indices.

Return type:

tuple[ndarray, ndarray, tuple[ndarray, ndarray]]

Returns:

Tuple containing: Tuple containing: fractional coordinates for centers and framework atoms

and Tuple containing: indices for centers used in the calculation of the diffusion and indices of framework atoms.

kinisi.parser.is_subset_approx(B, A, tol=1e-09)#

Check if all elements in B are approximately equal to any element in A within a tolerance. This is useful for comparing floating-point numbers where exact equality is not feasible.

Parameters:
  • B (array) – The array to check if it is a subset of A.

  • A (array) – The array to check against.

  • tol (float) – The tolerance for comparison. Default is 1e-9.

Return type:

bool

Returns:

True if all elements in B are approximately equal to any element in A, False otherwise.

kinisi.parser.is_orthorhombic(latt)#

Check if trajectory is always orthorhombic.

This function works by flattening each frames lattice vectors, then checking which are close to 0, and counting how many return True. If the cell is orthorhombic, only 3 elements, of 9, should be nonzero, leaving 6. Hence, count_nonzero should equal 6 on every element to return true. This does not measure lattice angles and so all vectors must be aligned with axes to return true.

Parameters:

latt (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – a scipp array with dimensions time,`dimension1`, and dimension2.

Return type:

bool

Returns:

True if lattice vectors are orthorhombic for all trajectory frames.

The ASEParser class is a parser for ASE atoms object. It is used to extract the necessary data for diffusion analysis from ASE.

class kinisi.ase.ASEParser(atoms, specie, time_step, step_skip, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, drift_indices=None, masses=None, progress=True)#

Bases: Parser

A parser for ASE Atoms objects

Parameters:
  • atoms (list[ase.atoms.Atoms]) – Atoms ordered in sequence of run.

  • specie (Union[str, list, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]) – symbol to calculate diffusivity for as a String, list of strings, or scipp.Variable of strings.

  • time_step (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Time step, in picoseconds, between steps in trajectory.

  • step_skip (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Sampling frequency of the trajectory (time_step is multiplied by this number to get the real time between output from the simulation file).

  • sub_sample_traj – Multiple of the time_step to sub sample at. Optional, defaults to 1 where all timesteps are used.

  • min_dt – Minimum time interval to be evaluated, in the simulation units. Optional, defaults to the produce of time_step and step_skip.

  • max_dt – Maximum time interval to be evaluated, in the simulation units. Optional, defaults to the length of the simulation.

  • n_steps – Number of steps to be used in the time interval function. Optional, defaults to 100 unless this is fewer than the total number of steps in the trajectory when it defaults to this number.

  • spacing – The spacing of the steps that define the time interval, can be either 'linear' or 'logarithmic'. If 'logarithmic' the number of steps will be less than or equal to that in the n_steps argument, such that all values are unique. Optional, defaults to linear.

  • sampling – The ways that the time-windows are sampled. The options are 'single-origin' or 'multi-origin' with the former resulting in only one observation per atom per time-window and the latter giving the maximum number of origins without sampling overlapping trajectories. Optional, defaults to 'multi-origin'.

  • memory_limit – Upper limit in the amount of computer memory that the displacements can occupy in gigabytes (GB). Optional, defaults to 8..

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

  • specie_indices (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Optional, list of indices to calculate diffusivity for as a list of indices. Specie must be set to None for this to function. Molecules can be specificed as a list of lists of indices. The inner lists must all be on the same length.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Optional, list of masses associated with the indices in specie_indices. Must be same shape as specie_indices.

  • framework_indices – Optional, list of framework indices to be used to correct framework drift. If an empty list is passed no drift correction will be performed.

get_structure_coords_latt(atoms, distance_unit, progress=True)#

Obtain the initial structure and displacement from a list of :py:class`pymatgen.core.structure.Structure`.

Parameters:
  • structures – Structures ordered in sequence of run.

  • sub_sample_traj – Multiple of the time_step to sub sample at. Optional, default is 1.

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

Return type:

tuple[ase.atoms.Atoms, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]), TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]

Returns:

Tuple containing: initial structure, fractional coordinates for all atoms, and lattice descriptions.

get_indices(structure, specie)#

Determine framework and non-framework indices for a ase compatible file.

Parameters:
  • structure (ase.atoms.Atoms) – Initial structure.

  • specie (Union[str, list, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]) – Specie to calculate diffusivity for as a String, list of strings, or scipp.Variable of strings.

  • framework_indices – Indices of framework to be used in drift correction. If set to None will return all indices that are not in indices.

Return type:

tuple[TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]), TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]

Returns:

Tuple containing: indices for the atoms in the trajectory used in the calculation of the diffusion and indices of framework atoms.

get_drift_indices(structure, specie_indices)#

Determine framework indices for an ase structure.

Parameters:

structure (ase.atoms.Atoms) – Initial structure.

Return type:

TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])

Returns:

Indices for the atoms in the trajectory used as framework atoms.

The MDAnalysisParser class is a parser for MDAnalysis universe object. It is used to extract the necessary data for diffusion analysis from an MDAnalysis universe.

class kinisi.mdanalysis.MDAnalysisParser(universe, specie, time_step, step_skip, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, drift_indices=None, masses=None, progress=True)#

Bases: Parser

Parser for MDAnalysis structures.

Takes an MDAnalysis.Universe object as an input.

Parameters:
  • universe (MDAnalysis.core.universe.Universe) – MDanalysis universe object to be parsed

  • specie (Union[str, list, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]) – Specie to calculate diffusivity for as a String, list of strings, or scipp.Variable of strings.

  • 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.

  • 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 used in the input structures. Optional, defaults to angstroms.

  • sub_sample_atoms – Subsample the atoms in the trajectory. Optional, defaults to 1.

  • sub_sample_traj – Subsample the trajectory. Optional, defaults to 1.

  • progress (bool) – Whether to show a progress bar when reading in the structures. Optional, defaults to True.

get_structure_coords_latt(universe, distance_unit, progress=True)#

Obtain the initial structure, coordinates, and lattice parameters from an MDAnalysis.Universe object.

Parameters:
  • universe (MDAnalysis.core.universe.Universe) – MDanalysis universe object.

  • progress (bool) – Whether to show a progress bar when reading in the structures.

  • sub_sample_atoms – Subsample the atoms in the trajectory. Optional, defaults to 1.

  • sub_sample_traj – Subsample the trajectory. Optional, defaults to 1.

Return type:

tuple[MDAnalysis.core.universe.Universe, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]), TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]

Returns:

A tuple of: the initial structure (as a MDAnalysis.core.universe.Universe.atoms), coordinates (as a scipp array with dimensions of time, atom, and dimension), and lattice parameters (as a scipp array with dimensions time, dimension1, and dimension2).

get_indices(structure, specie)#

Determine framework and non-framework indices for an MDAnalysis compatible file.

Parameters:
  • structure (MDAnalysis.universe.Universe.atoms) – Initial structure.

  • specie (Union[str, list, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]) – Specie to calculate diffusivity for as a String, list of strings, or scipp.Variable of strings.

Return type:

tuple[TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]), TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]

Returns:

Tuple containing indices for the atoms in the trajectory used in the calculation of the diffusion and indices of framework atoms.

get_drift_indices(structure, specie_indices)#

Determine framework indices for an MDAnalysis structure.

Parameters:

structure (MDAnalysis.universe.Universe.atoms) – Initial structure.

Return type:

TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])

Returns:

Indices for the atoms in the trajectory used as framework atoms.

The PymatgenParser class is a parser for pymatgen structures for kinisi. It is used to extract the necessary data for diffusion analysis from a list of pymatgen structures.

class kinisi.pymatgen.PymatgenParser(structures, specie, time_step, step_skip, dt=None, dimension='xyz', distance_unit=Unit(1e-10 * m), specie_indices=None, drift_indices=None, masses=None, progress=True)#

Bases: Parser

Parser for pymatgen structures.

This takes a list of pymatgen structures as an input.

Parameters:
  • structures (list[pymatgen.core.structure.Structure]) – Structures ordered in sequence of run.

  • specie (Union[str, list, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]) – Specie to calculate diffusivity for as a String, list of strings, or scipp.Variable of strings.

  • 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.

  • 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 used in the input structures. Optional, defaults to angstroms.

  • specie_indices (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Indices of the specie to calculate the diffusivity for. Optional, defaults to None.

  • masses (TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])) – Masses of the atoms in the structure. Optional, defaults to None.

  • progress (bool) – Whether to show a progress bar when reading in the structures. Optional, defaults to True.

get_structure_coords_latt(structures, distance_unit, progress=True)#

Obtain the initial structure, coordinates, and lattice parameters from a list of pymatgen structures.

Parameters:
  • structures (list[pymatgen.core.structure.Structure]) – Structures ordered in sequence of run.

  • progress (bool) – Whether to show a progress bar when reading in the structures.

Return type:

tuple[pymatgen.core.structure.Structure, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]), TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]

Returns:

A tuple of the initial structure (as a pymatgen.core.structure.Structure), coordinates (as a scipp array with dimensions of time, atom, and dimension), and lattice parameters (as a scipp array with dimensions time, dimension1, and dimension2).

get_indices(structure, specie)#

Determine the framework and mobile indices from a pymatgen structure.

Parameters:
  • structure (pymatgen.core.structure.Structure) – The initial structure to determine the indices from.

  • specie (Union[str, list, TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]) – The specie to calculate the diffusivity for as a String, list of strings, or scipp.Variable of strings.

Return type:

tuple[TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any]), TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])]

Returns:

A tuple of the indices for the specie of interest (mobile) and the drift (framework) indices.

get_drift_indices(structure, specie_indices)#

Determine framework indices for an pymatgen structure.

Parameters:

structure (pymatgen.core.structure.Structure) – Initial structure.

Return type:

TypeVar(VariableLikeType, Variable, DataArray, Dataset, DataGroup[Any])

Returns:

Indices for the atoms in the trajectory used as framework atoms.