mdgo.conductivity module

This module implements functions to calculate the ionic conductivity.

mdgo.conductivity.autocorr_fft(x)[source]

Calculates the autocorrelation function using the fast Fourier transform.

Parameters

x (numpy.array) – function on which to compute autocorrelation function

Return type

numpy.ndarray

Returns a numpy.array of the autocorrelation function

Return type

ndarray

Parameters

x (numpy.ndarray) –

mdgo.conductivity.msd_fft(r)[source]

Calculates mean square displacement of the array r using the fast Fourier transform.

Parameters

r (numpy.array) – atom positions over time

Return type

numpy.ndarray

Returns a numpy.array containing the mean-squared displacement over time

Return type

ndarray

Parameters

r (numpy.ndarray) –

mdgo.conductivity.calc_cond_msd(u, anions, cations, run_start, cation_charge=1, anion_charge=-1)[source]

Calculates the conductivity “mean square displacement” over time

Note

Coordinates must be unwrapped (in dcd file when creating MDAnalysis Universe) Ions selections may consist of only one atom per ion, or include all of the atoms

in the ion. The ion AtomGroups may consist of multiple types of cations/anions.

Parameters
  • u (Universe) – MDAnalysis universe

  • anions (AtomGroup) – MDAnalysis AtomGroup containing all anions

  • cations (AtomGroup) – MDAnalysis AtomGroup containing all cations

  • run_start (int) – index of trajectory from which to start analysis

  • cation_charge (int) – net charge of cation

  • anion_charge (int) – net charge of anion

Return type

numpy.ndarray

Returns a numpy.array containing conductivity “MSD” over time

Return type

ndarray

Parameters
  • u (MDAnalysis.Universe) –

  • anions (MDAnalysis.AtomGroup) –

  • cations (MDAnalysis.AtomGroup) –

  • run_start (int) –

  • cation_charge (Union[int, float]) –

  • anion_charge (Union[int, float]) –

mdgo.conductivity.get_beta(msd, time_array, start, end)[source]

Fits the MSD to the form t^(beta) and returns beta. beta = 1 corresponds to the diffusive regime.

Parameters
  • msd (numpy.array) – mean squared displacement

  • time_array (numpy.array) – times at which position data was collected in the simulation

  • start (int) – index at which to start fitting linear regime of the MSD

  • end (int) – index at which to end fitting linear regime of the MSD

Return type

tuple

Returns beta (int) and the range of beta values within the region

Return type

tuple

Parameters
  • msd (numpy.ndarray) –

  • time_array (numpy.ndarray) –

  • start (int) –

  • end (int) –

mdgo.conductivity.choose_msd_fitting_region(msd, time_array)[source]

Chooses the optimal fitting regime for a mean-squared displacement. The MSD should be of the form t^(beta), where beta = 1 corresponds to the diffusive regime; as a rule of thumb, the MSD should exhibit this linear behavior for at least a decade of time. Finds the region of the MSD with the beta value closest to 1.

Note

If a beta value great than 0.9 cannot be found, returns a warning that the computed conductivity may not be reliable, and that longer simulations or more replicates are necessary.

Parameters
  • msd (numpy.array) – mean squared displacement

  • time_array (numpy.array) – times at which position data was collected in the simulation

Return type

tuple

Returns at tuple with the start of the fitting regime (int), end of the fitting regime (int), and the beta value of the fitting regime (float).

Return type

tuple

Parameters
  • msd (numpy.ndarray) –

  • time_array (numpy.ndarray) –

mdgo.conductivity.conductivity_calculator(time_array, cond_array, v, name, start, end, T, units='real')[source]

Calculates the overall conductivity of the system

Parameters
  • time_array (numpy.array) – times at which position data was collected in the simulation

  • cond_array (numpy.array) – conductivity “mean squared displacement”

  • v (float) – simulation volume (Angstroms^3)

  • name (str) – system name

  • start (int) – index at which to start fitting linear regime of the MSD

  • end (int) – index at which to end fitting linear regime of the MSD

  • units (str) – unit system (currently ‘real’ and ‘lj’ are supported)

  • T (Union[int, float]) –

Return type

float

Returns the overall ionic conductivity (float)

Return type

float

Parameters
  • time_array (numpy.ndarray) –

  • cond_array (numpy.ndarray) –

  • v (Union[int, float]) –

  • name (str) –

  • start (int) –

  • end (int) –

  • T (Union[int, float]) –

  • units (str) –