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.DiffusionAnalyzer class enables the analysis of mean-squared displacement and the self-diffusion coefficient.

class kinisi.diffusion_analyzer.DiffusionAnalyzer(delta_t, disp_3d, n_o, volume)#

Bases: Analyzer

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

Parameters:
  • delta_t (ndarray) – An array of the timestep values.

  • disp_3d (List[ndarray]) – A list of arrays, where each array has the axes [atom, displacement observation, dimension]. There is one array in the list for each delta_t value. Note: it is necessary to use a list of arrays as the number of observations is not necessary the same at each data point.

  • volume (float) – The volume of the simulation cell.

  • uncertainty_params – The parameters for the kinisi.diffusion.DiffBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

to_dict()#
Return type:

dict

Returns:

Dictionary description of DiffusionAnalyzer.

classmethod from_dict(my_dict)#

Generate a DiffusionAnalyzer object from a dictionary.

Parameters:

my_dict (dict) – The input dictionary.

Return type:

DiffusionAnalyzer

Returns:

New DiffusionAnalyzer object.

classmethod from_pymatgen(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a DiffusionAnalyzer object from a list or nested list of pymatgen.core.structure.Structure objects.

Parameters:
  • trajectory (List[Union[pymatgen.core.structure.Structure, List[pymatgen.core.structure.Structure]]]) – The list or nested list of structures to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – If trajectory is a list of pymatgen.core.structure.Structure objects, this should be None. However, if a list of lists 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Dtype trajectory:

List[Union[‘pymatgen.core.structure.Structure’, List[‘pymatgen.core.structure.Structure’]]]

Returns:

Relevant DiffusionAnalyzer object.

classmethod from_ase(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a DiffusionAnalyzer object from a list or nested list of ase.atoms.Atoms objects.

Parameters:
  • trajectory (List[Union[ase.atoms.Atoms, List[ase.atoms.Atoms]]]) – The list or nested list of structures to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.ASEParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – If trajectory is a list of ase.atoms.Atoms objects, this should be None. However, if a list of lists 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Dtype trajectory:

List[Union[‘ase.atoms.Atoms’, List[‘ase.atoms.Atoms’]]]

Returns:

Relevant DiffusionAnalyzer object.

classmethod from_Xdatcar(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a DiffusionAnalyzer object 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 Xdatcar or list of Xdatcar objects to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant DiffusionAnalyzer object.

classmethod from_file(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a DiffusionAnalyzer object from a single or a list of Xdatcar file(s).

Parameters:
  • trajectory (Union[str, List[str]]) – The file or list of Xdatcar files to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (Optional[str]) – If trajectory is a single file, this should be None. However, if a list of files 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.

  • uncertainty_params (Optional[dict]) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant DiffusionAnalyzer object.

classmethod from_universe(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create an DiffusionAnalyzer object from an MDAnalysis.core.universe.Universe object.

Parameters:
  • trajectory (MDAnalysis.core.universe.Universe) – The universe to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.MDAnalysisParser object. See the appropriate documention for more guidance on this dictionary.

  • dtype (str) – If trajectory is a single file, this should be None. However, if a list of files is passed, then it is necessary to identify that these constitute a series of identical starting points with different random seeds, in which case the dtype should be identical. For a series of consecutive trajectories, please construct the relevant object using MDAnalysis.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant DiffusionAnalyzer object.

diffusion(start_dt, diffusion_params=None)#

Calculate the diffusion coefficicent using the bootstrap-GLS methodology.

Parameters:
  • start_dt (float) – The starting time for the analysis to find the diffusion coefficient. This should be the start of the diffusive regime in the simulation.

  • diffusion_params (Optional[dict]) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

property msd: ndarray#
Returns:

MSD for the input trajectories. Note that this is the bootstrap sampled MSD, not the numerical

average from the data.

property msd_std: ndarray#
Returns:

MSD standard deviation values for the input trajectories (a single standard deviation).

property D: uravu.distribution.Distribution#
Returns:

Diffusion coefficient distribution.

property flatchain: ndarray#
Returns:

sampling flatchain

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

class kinisi.jump_diffusion_analyzer.JumpDiffusionAnalyzer(delta_t, disp_3d, n_o, volume)#

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:
  • delta_t (ndarray) – An array of the timestep values.

  • disp_3d (List[ndarray]) – A list of arrays, where each array has the axes [atom, displacement observation, dimension]. There is one array in the list for each delta_t value. Note: it is necessary to use a list of arrays as the number of observations is not necessary the same at each data point.

  • volume (float) – The volume of the simulation cell.

  • uncertainty_params – The parameters for the kinisi.diffusion.DiffBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

to_dict()#
Return type:

dict

Returns:

Dictionary description of JumpDiffusionAnalyzer.

classmethod from_dict(my_dict)#

Generate a DiffusionAnalyzer object from a dictionary.

Parameters:

my_dict – The input dictionary.

Return type:

JumpDiffusionAnalyzer

Returns:

New DiffusionAnalyzer object.

classmethod from_pymatgen(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a JumpDiffusionAnalyzer object from a list or nested list of pymatgen.core.structure.Structure objects.

Parameters:
  • trajectory (List[Union[pymatgen.core.structure.Structure, List[pymatgen.core.structure.Structure]]]) – The list or nested list of structures to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – If trajectory is a list of pymatgen.core.structure.Structure objects, this should be None. However, if a list of lists 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant JumpDiffusionAnalyzer object.

classmethod from_ase(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a JumpDiffusionAnalyzer object from a list or nested list of ase.atoms.Atoms objects.

Parameters:
  • trajectory (List[Union[ase.atoms.Atoms, List[ase.atoms.Atoms]]]) – The list or nested list of structures to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.AseParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – If trajectory is a list of ase.atoms.Atoms objects, this should be None. However, if a list of lists 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant JumpDiffusionAnalyzer object.

classmethod from_Xdatcar(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a JumpDiffusionAnalyzer object 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 Xdatcar or list of Xdatcar objects to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant JumpDiffusionAnalyzer object.

classmethod from_file(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create a JumpDiffusionAnalyzer object from a single or a list of Xdatcar file(s).

Parameters:
  • trajectory (Union[str, List[str]]) – The file or list of Xdatcar files to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (Optional[str]) – If trajectory is a single file, this should be None. However, if a list of files 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.

  • uncertainty_params (Optional[dict]) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant JumpDiffusionAnalyzer object.

classmethod from_universe(trajectory, parser_params, dtype=None, uncertainty_params=None)#

Create an JumpDiffusionAnalyzer object from an MDAnalysis.core.universe.Universe object.

Parameters:
  • trajectory (MDAnalysis.core.universe.Universe) – The universe to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.MDAnalysisParser object. See the appropriate documention for more guidance on this dictionary.

  • dtype (str) – If trajectory is a single file, this should be None. However, if a list of files is passed, then it is necessary to identify that these constitute a series of identical starting points with different random seeds, in which case the dtype should be identical. For a series of consecutive trajectories, please construct the relevant object using MDAnalysis.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

Returns:

Relevant JumpDiffusionAnalyzer object.

jump_diffusion(start_dt, jump_diffusion_params=None)#

Calculate the jump diffusion coefficicent using the bootstrap-GLS methodology.

Parameters:
  • start_dt (float) – The starting time for the analysis to find the diffusion coefficient. This should be the start of the diffusive regime in the simulation.

  • ump_diffusion_params – The parameters for the kinisi.diffusion.MSTDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

property mstd: ndarray#
Returns:

MSTD for the input trajectories. Note that this is the bootstrap sampled MSD, not the numerical average from the data.

property mstd_std: ndarray#
Returns:

MSD standard deviation values for the input trajectories (a single standard deviation).

property D_J: uravu.distribution.Distribution#
Returns:

Jump diffusion coefficient

property flatchain: ndarray#
Returns:

sampling flatchain

The kinisi.analyze.ConductivityAnalyzer class will allow the conductivity of a material in a simulation to be found, without assuming a Haven ration of 1.

class kinisi.conductivity_analyzer.ConductivityAnalyzer(delta_t, disp_3d, n_o, volume)#

Bases: Analyzer

The kinisi.analyze.ConductivityAnalyzer class performs analysis of conductive relationships in materials. This is achieved through the application of a bootstrapping methodology to obtain the most statistically accurate values for mean squared charge displacement uncertainty and covariance. The time-dependence of the MSCD 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:
  • delta_t (ndarray) – An array of the timestep values.

  • disp_3d (List[ndarray]) – A list of arrays, where each array has the axes [atom, displacement observation, dimension]. There is one array in the list for each delta_t value. Note: it is necessary to use a list of arrays as the number of observations is not necessary the same at each data point.

  • volume (float) – The volume of the simulation cell.

  • uncertainty_params – The parameters for the kinisi.diffusion.DiffBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

  • ionic_charge – The charge on the mobile ions, either an array with a value for each ion or a scalar if all values are the same. Optional, default is 1.

to_dict()#
Return type:

dict

Returns:

Dictionary description of ConductivityAnalyzer.

classmethod from_dict(my_dict)#

Generate a ConductivityAnalyzer object from a dictionary.

Parameters:

my_dict (dict) – The input dictionary.

Return type:

ConductivityAnalyzer

Returns:

New :py:class`ConductivityAnalyzer` object.

classmethod from_pymatgen(trajectory, parser_params, dtype=None, uncertainty_params=None, ionic_charge=1)#

Create a ConductivityAnalyzer object from a list or nested list of pymatgen.core.structure.Structure objects.

Parameters:
  • trajectory (List[Union[pymatgen.core.structure.Structure, List[pymatgen.core.structure.Structure]]]) – The list or nested list of structures to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – If trajectory is a list of pymatgen.core.structure.Structure objects, this should be None. However, if a list of lists 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

  • ionic_charge (Union[ndarray, int]) – The charge on the mobile ions, either an array with a value for each ion or a scalar if all values are the same. Optional, default is 1.

Returns:

Relevant ConductivityAnalyzer object.

classmethod from_ase(trajectory, parser_params, dtype=None, uncertainty_params=None, ionic_charge=1)#

Create a ConductivityAnalyzer object from a list or nested list of ase.Atoms objects.

Parameters:
  • trajectory (List[Union[ase.atoms.Atoms, List[ase.atom.Atoms]]]) – The list or nested list of structures to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.ASEParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – If trajectory is a list of ase.Atoms objects, this should be None. However, if a list of lists 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

  • ionic_charge (Union[ndarray, int]) – The charge on the mobile ions, either an array with a value for each ion or a scalar if all values are the same. Optional, default is 1.

Returns:

Relevant ConductivityAnalyzer object.

classmethod from_Xdatcar(trajectory, parser_params, dtype=None, uncertainty_params=None, ionic_charge=1)#

Create a ConductivityAnalyzer object 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 Xdatcar or list of Xdatcar objects to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (str) – 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.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

  • ionic_charge (Union[ndarray, int]) – The charge on the mobile ions, either an array with a value for each ion or a scalar if all values are the same. Optional, default is 1.

Returns:

Relevant ConductivityAnalyzer object.

classmethod from_file(trajectory, parser_params, dtype=None, uncertainty_params=None, ionic_charge=1)#

Create a ConductivityAnalyzer object from a single or a list of Xdatcar file(s).

Parameters:
  • trajectory (Union[str, List[str]]) – The file or list of Xdatcar files to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.PymatgenParser object. See the appropriate documentation for more guidance on this dictionary.

  • dtype (Optional[str]) – If trajectory is a single file, this should be None. However, if a list of files 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.

  • uncertainty_params (Optional[dict]) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

  • ionic_charge (Union[ndarray, int]) – The charge on the mobile ions, either an array with a value for each ion or a scalar if all values are the same. Optional, default is 1.

Returns:

Relevant ConductivityAnalyzer object.

classmethod from_universe(trajectory, parser_params, dtype=None, uncertainty_params=None, ionic_charge=1)#

Create an ConductivityAnalyzer object from an MDAnalysis.core.universe.Universe object.

Parameters:
  • trajectory (MDAnalysis.core.universe.Universe) – The universe to be analysed.

  • parser_params (dict) – The parameters for the kinisi.parser.MDAnalysisParser object. See the appropriate documention for more guidance on this dictionary.

  • dtype (str) – If trajectory is a single file, this should be None. However, if a list of files is passed, then it is necessary to identify that these constitute a series of identical starting points with different random seeds, in which case the dtype should be identical. For a series of consecutive trajectories, please construct the relevant object using MDAnalysis.

  • uncertainty_params (dict) – The parameters for the kinisi.diffusion.MSDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters.

  • ionic_charge (Union[ndarray, int]) – The charge on the mobile ions, either an array with a value for each ion or a scalar if all values are the same. Optional, default is 1.

Returns:

Relevant ConductivityAnalyzer object.

conductivity(start_dt, temperature, conductivity_params=None)#

Calculate the jump diffusion coefficicent using the bootstrap-GLS methodology.

Parameters:
  • start_dt (float) – The starting time for the analysis to find the diffusion coefficient. This should be the start of the diffusive regime in the simulation.

  • temperature (float) – Simulation temperature in Kelvin

  • conductivity_params (Optional[dict]) – The parameters for the kinisi.diffusion.MSCDBootstrap object. See the appropriate documentation for more guidance on this. Optional, default is the default bootstrap parameters

property mscd: ndarray#
Returns:

MSCD for the input trajectories. Note that this is the bootstrap sampled value, not the numerical average from the data.

property mscd_std: ndarray#
Returns:

MSCD standard deviation values for the input trajectories (a single standard deviation).

property sigma: uravu.distribution.Distribution#
Returns:

Conductivity, in mS^{1}cm^{-1}.

property flatchain: ndarray#
Returns:

sampling flatchain

This module contains the base class for the different Analyzer objects used by kinisi.

class kinisi.analyzer.Analyzer(delta_t, disp_3d, n_o, volume)#

Bases: object

This class is the superclass for the kinisi.analyze.DiffusionAnalyzer, kinisi.analyze.JumpDiffusionAnalyzer and kinisi.analyze.ConductivityAnalyzer classes. Therefore all of the properties here are available to these other classes.

Parameters:
  • delta_t (ndarray) – An array of the timestep values.

  • disp_3d (List[ndarray]) – A list of arrays, where each array has the axes [atom, displacement observation, dimension]. There is one array in the list for each delta_t value. Note: it is necessary to use a list of arrays as the number of observations is not necessary the same at each data point.

  • volume (float) – The volume of the simulation cell.

save(filename)#

Save the Analyzer object as a HDF5 file.

Parameters:

filename (str) – Name for the file, no file extension is required and if one if given it is replaced with .hdf.

classmethod load(filename)#

Load the Analyzer object from an HDF5 file.

Parameters:

filename (str) – Name for the file, any file extension will be replaced with .hdf.

Return type:

Analyzer

Returns:

An Analyzer object from the file.

to_dict()#
Return type:

dict

Returns:

Dictionary description of Analyzer.

classmethod from_dict(my_dict)#

Generate an Analyzer object from a dictionary.

Parameters:

my_dict (dict) – The input dictionary.

Return type:

Analyzer

Returns:

New Analyzer object.

posterior_predictive(posterior_predictive_params=None)#

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 distribution: ndarray#
Returns:

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

plotting of credible intervals.

property dt: ndarray#
Returns:

Timestep values that have been sampled.

property dr: List[uravu.distribution.Distribution]#
Returns:

A list of uravu.distribution.Distribution objects that describe the Euclidean displacement at each dt.

property ngp_max: float#
Returns:

Position in dt where the non-Gaussian parameter is maximised.

property intercept: uravu.distribution.Distribution#
Returns:

The distribution describing the intercept.

property volume: float#
Returns:

Volume of system, in cubic angstrom.