Onsager

class mdhelper.analysis.transport.Onsager(groups: AtomGroup | tuple[AtomGroup], groupings: str | tuple[str] = 'atoms', temperature: float | Quantity | Quantity = 300, *, charges: ndarray[float] | Quantity | Quantity = None, dimensions: ndarray[float] | Quantity | Quantity = None, dt: float | Quantity | Quantity = None, n_blocks: int = 1, center: bool = False, center_atom: bool = False, center_wrap: bool = False, fft: bool = True, reduced: bool = False, unwrap: bool = False, verbose: bool = True, **kwargs)[source]

Bases: SerialAnalysisBase

A serial implementation to calculate the Onsager transport coefficients and its related properties.

Note

The simulation must have been run with a constant timestep \(\Delta t\) and the frames to be analyzed must be evenly spaced and proceed forward in time for this analysis module to function correctly.

The Onsager transport framework [1] can be used to analyze transport properties in bulk constant-volume fluids and electrolytic systems, such as those in the canonical (\(NVT\)), microcanonical (\(NVE\)), and grand canonical (\(\mu VT\)) ensembles.

The Onsager transport equation

\[\pmb{J}_i=-\sum_j L_{ij}\nabla\bar{\mu}_j\]

relates the flux \(\pmb{J}_i\) of species \(i\) to the Onsager transport coefficients \(L_{ij}\) and the electrochemical potential \(\bar{\mu}_j\) of species \(j\). There is an Onsager transport coefficient for each pair of species that, unlike the Nernst–Einstein equation, captures the strong cross-correlations in electrolytes.

The Onsager transport coefficients can be calculated from the particle positions over time using the Einstein relation

\[L_{ij}=\frac{1}{6k_\mathrm{B}TV} \lim_{t\rightarrow\infty}\frac{d}{dt}\left\langle\sum_\alpha \left[\pmb{r}_{i,\alpha}(t)-\pmb{r}_{i,\alpha}(0)\right] \cdot\sum_\beta\left[\pmb{r}_{j,\beta}(t) -\pmb{r}_{j,\beta}(0)\right]\right\rangle\]

where \(k_\mathrm{B}\) is the Boltzmann constant, \(T\) is the system temperature, \(V\) is the system volume, \(t\) is time, and \(\pmb{r}_\alpha\) and \(\pmb{r}_\beta\) are the positions of particles \(\alpha\) and \(\beta\) belonging to species \(i\) and \(j\), respectively. The angular brackets denote the ensemble average. It is evident that \(L_{ij}=L_{ji}\); hence, the equation above is an Onsager reciprocal relation.

The diagonal terms \(L_{ii}\) in the matrix of Onsager transport coefficients captures the self-diffusion and self-correlations for a single species \(i\) and has two contributions:

\[L_{ii}=L_{ii}^\mathrm{self}+L_{ii}^\mathrm{distinct}\]

The self term

\[L_{ii}^\mathrm{self}=\frac{1}{6k_\mathrm{B}TV} \lim_{t\rightarrow\infty}\frac{d}{dt} \sum_\alpha\left\langle\left[ \pmb{r}_{i,\alpha}(t)-\pmb{r}_{i,\alpha}(0) \right]^2\right\rangle\]

is given by the autocorrelation function of the flux of a single particle \(\alpha\) when \(\alpha=\beta\), while the distinct term

\[L_{ii}^\mathrm{self}=\frac{1}{6k_\mathrm{B}TV} \lim_{t\rightarrow\infty}\frac{d}{dt} \sum_\alpha\sum_{\beta\neq\alpha}\left\langle\left[ \pmb{r}_{i,\alpha}(t)-\pmb{r}_{i,\alpha}(0)\right]\cdot \left[\pmb{r}_{i,\beta}(t)-\pmb{r}_{i,\beta}(0)\right] \right\rangle\]

is given by the cross-correlation between two distinct particles when \(\alpha\neq\beta\). Naturally, the self term is related to the self-diffusion coefficient \(D_i\) via

\[L_{ii}^\mathrm{self}=\frac{D_i\rho_i}{k_\mathrm{B}T}\]

where \(\rho_i\) is the number density of species \(i\).

Finally, the Onsager transport coefficients can be used to obtain experimentally relevant transport properties:

  • Ionic conductivity:

    \[\kappa=F^2\sum_i\sum_j z_iz_jL_{ij}\]

    \(F\) is the Faraday constant.

  • Electrophoretic mobility:

    \[\mu_i=\frac{F}{\rho_i}\sum_j z_jL_{ij}\]
  • Transference number:

    \[t_i=\frac{\sum_j z_iz_jL_{ij}}{\sum_k\sum_l z_kz_lL_{kl}}\]
Parameters:
groupsMDAnalysis.AtomGroup or array-like

Group(s) of atoms for which the mean squared displacements (or analogs) are calculated.

groupingsstr or array-like, default: "atoms"

Determines whether the centers of mass are used in lieu of individual atom positions. If groupings is a str, the same value is used for all groups.

Note

In a standard trajectory file, segments (or chains) contain residues (or molecules), and residues contain atoms. This heirarchy must be adhered to for this analysis module to function correctly, unless your selected grouping is always "atoms".

Valid values:

  • "atoms": Atom positions (generally for coarse-grained simulations).

  • "residues": Residues’ centers of mass (for atomistic simulations).

  • "segments": Segments’ centers of mass (for atomistic polymer simulations).

temperaturefloat, openmm.unit.Quantity, or pint.Quantity, default: 300

System temperature \(T\).

Note

If reduced=True, temperature should be equal to the energy scale. When the Lennard-Jones potential is used, it generally means that \(T^* = 1\), or temperature=1.

Reference unit: \(\mathrm{K}\).

chargesarray-like, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional

Charge numbers \(z_i\) for the specified groupings in the \(N_\mathrm{g}\) groups. If not provided, it will be retrieved from the topology or trajectory.

Shape: \((N_\mathrm{g},)\).

Reference unit: \(\mathrm{e}\).

dimensionsarray-like, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional

System dimensions. If not provided, they are retrieved from the topology or trajectory.

Tip

You can zero out dimensions by setting them to 0 if your simulation is pseudo-1D or pseudo-2D. For example, if you have a one-layer thick slab in the \(xy\)-plane, you can set dimensions=np.array((L_x, L_y, 0)) to evaluate the transport properties in 2D. Note that this affects the reference units of the Onsager transport coefficients and its related properties, but not the self-diffusivity.

Shape: \((3,)\).

Reference unit: \(\mathrm{Å}\).

dtfloat, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional

Time between frames \(\Delta t\). While this is normally determined from the trajectory, the trajectory may not have the correct information if the data is in reduced units. For example, if your reduced timestep is \(0.01\) and you output trajectory data every \(10,000\) timesteps, then \(\Delta t=100\).

Reference unit: \(\mathrm{ps}\).

n_blocksint, keyword-only, default: 1

Number of blocks \(N_\mathrm{b}\) to split the trajectory into.

centerbool, keyword-only, default: False

Determines whether the system center of mass is subtracted from the positions used to compute the transport properties.

center_atombool, keyword-only, default: False

Determines whether the system center of mass is computed using the atom positions, regardless of groupings. Has no effect if center=False.

center_wrapbool, keyword-only, default: False

Determines whether the system center of mass is computed using the wrapped particle positions. Has no effect if center=False.

fftbool, keyword-only, default: True

Determines whether fast Fourier transforms (FFT) are used to evaluate the mean squared displacements (or analogs).

reducedbool, keyword-only, default: False

Specifies whether the data is in reduced units. Affects temperature, results.times, etc.

unwrapbool, keyword-only, default: False

Determines if atom positions are unwrapped. Ensure that unwrap=False when the trajectory already contains unwrapped particle positions, as this parameter is used in conjunction with center_wrap to determine the appropriate system center of mass.

verbosebool, keyword-only, default: True

Determines whether detailed progress is shown.

**kwargs

Additional keyword arguments to pass to MDAnalysis.analysis.base.AnalysisBase.

Notes

  • The values in results.msd_cross are actually “summed squared displacements”—that is, results.msd_cross contains the sum of all particles’ squared displacements instead of their average. However, it is extremely unlikely the raw values in results.msd_cross will be used other than to calculate the Onsager transport coefficients. On the other hand, results.msd_self is averaged over all particles. Therefore, it can be plotted against results.times to get the self-diffusion coefficients.

References

[1]

Fong, K. D.; Self, J.; McCloskey, B. D.; Persson, K. A. Onsager Transport Coefficients and Transference Numbers in Polyelectrolyte Solutions and Polymerized Ionic Liquids. Macromolecules 2020, 53 (21), 9503–9512. https://doi.org/10.1021/acs.macromol.0c02001.

Attributes:
universeMDAnalysis.Universe

MDAnalysis.core.universe.Universe object containing all information describing the system.

results.unitsdict

Reference units for the results. For example, to get the reference units for results.times, call results.units["results.times"].

results.pairstuple

All unique pairs of indices of the groups of atoms in groups. The ordering coincides with the column indices in results.msd.

results.timesnumpy.ndarray

Changes in time \(t-t_0\).

Shape: \((N_t,)\).

Reference unit: \(\mathrm{ps}\).

results.msd_crossnumpy.ndarray

MSDs (or analogs) for the \(N_\mathrm{g}\) groups over \(N_\mathrm{b}\) blocks of \(N_t\) trajectory frames each. Includes the dimensionality scaling factor. See Notes about what this value actually means.

Shape: \((_{N_\mathrm{g}+1}\mathrm{C}_{2},\,N_\mathrm{b},\, N_t)\).

Reference unit: \(\mathrm{Å}^2\).

results.msd_selfnumpy.ndarray

Self MSDs (or analogs) for the \(N_\mathrm{g}\) groups over \(N_\mathrm{b}\) blocks of \(N_t\) trajectory frames each. Includes the dimensionality scaling factor. See Notes about what this value actually means.

Shape: \((N_\mathrm{g},\,N_\mathrm{b},\,N_t)\).

Reference unit: \(\mathrm{Å}^2\).

results.L_ijnumpy.ndarray

Onsager transport coefficients \(L_{ij}\). Only available after running calculate_transport_coefficients().

Shape: \((N_\mathrm{b},\,N_\mathrm{g},\,N_\mathrm{g})\).

Reference unit: \(\mathrm{mol/(kJ}\cdot\mathrm{Å}\cdot\mathrm{ps)}\).

results.L_ii_selfnumpy.ndarray

Self Onsager transport coefficient terms \(L_{ii}^\mathrm{self}\). Note that \(L_{ii}^\mathrm{self}\) is related to \(D_i\) via

\[L_{ii}^\mathrm{self}=\dfrac{N}{k_\mathrm{B}TV}D_i\]

Only available after running calculate_transport_coefficients().

Shape: \((N_\mathrm{b},\,N_\mathrm{g})\).

Reference unit: \(\mathrm{mol/(kJ}\cdot\mathrm{Å}\cdot\mathrm{ps)}\).

results.D_inumpy.ndarray

Self-diffusion coefficients \(D_i\). Only available after running calculate_transport_coefficients().

Shape: \((N_\mathrm{b},\,N_\mathrm{g})\).

Reference unit: \(\mathrm{Å}^2/\mathrm{ps)}\).

results.conductivitiesnumpy.ndarray

Conductivities \(\kappa\). Only available after running calculate_conductivity().

Shape: \((N_\mathrm{b},\,N_\mathrm{g})\).

Reference unit: \(\mathrm{C}^2/(\mathrm{kJ}\cdot \mathrm{Å}\cdot\mathrm{ps})\).

To SI unit: \(1 imes10^{19}\,\mathrm{S}/\mathrm{m}\).

results.electrophoretic_mobilitiesnumpy.ndarray

Electrophoretic mobilities \(\mu_i\). Only available after running calculate_electrophoretic_mobility().

Shape: \((N_\mathrm{b},\,N_\mathrm{g})\).

Reference unit: \(\mathrm{Å}^2\cdot\mathrm{C/(kJ}\cdot\mathrm{ps)}\).

To SI unit: \(1 imes10^{-11}\,\mathrm{m}^2/ (\mathrm{V}\cdot\mathrm{s})\).

results.transference_numbersnumpy.ndarray

Transference numbers \(t_i\). Only available after running calculate_transference_number().

Shape: \((N_\mathrm{b},\,N_\mathrm{g})\).

Methods

calculate_conductivity

Calculates the ionic conductivity \(\kappa\) using the Onsager transport coefficients \(L_{ij}\).

calculate_electrophoretic_mobility

Calculates the electrophoretic mobility \(\mu_i\) of each species using the Onsager transport coefficients \(L_{ij}\).

calculate_transference_number

Calculates the transference number of each species using the Onsager transport coefficients \(L_{ij}\).

calculate_transport_coefficients

Fits the mean squared displacements (MSDs) (or analogs) to evaluate the Onsager transport coefficients \(L_{ij}\) and \(L_{ii}^\mathrm{self}\).

run

Performs the calculation in serial.

save

Saves results to a binary or archive file in NumPy format.

calculate_conductivity(*, charges: ndarray[float] | Quantity | Quantity = None) None[source]

Calculates the ionic conductivity \(\kappa\) using the Onsager transport coefficients \(L_{ij}\).

Parameters:
chargesarray-like, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional

Charge numbers \(z_i\) of the groupings in the \(N_\mathrm{g}\) groups. This argument is optional only if charges has previously been passed to a calculation method belonging to this class.

Shape: \((N_\mathrm{g},)\).

Reference unit: \(\mathrm{e}\).

calculate_electrophoretic_mobility(*, charges: ndarray[float] | Quantity | Quantity = None, rhos: ndarray[float] | Quantity | Quantity = None) None[source]

Calculates the electrophoretic mobility \(\mu_i\) of each species using the Onsager transport coefficients \(L_{ij}\).

Parameters:
chargesarray-like, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional

Charge numbers \(z_i\) of the groupings in the \(N_\mathrm{g}\) groups. This argument is optional only if charge information is present in the topology or charges has previously been passed to a calculation method belonging to this class.

Shape: \((N_\mathrm{g},)\).

Reference unit: \(\mathrm{e}\).

rhosarray-like, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional

Number densities \(n_i\) of the groupings in the \(N_\mathrm{g}\) groups.

Shape: \((N_\mathrm{g},)\).

Reference unit: \(\mathrm{Å}^{-3}\).

calculate_transference_number(*, charges: ndarray[float] | Quantity | Quantity = None) None[source]

Calculates the transference number of each species using the Onsager transport coefficients \(L_{ij}\).

Parameters:
chargesarray-like, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional

Charge numbers \(z_i\) of the groupings in the \(N_\mathrm{g}\) groups. This argument is optional only if charge information is present in the topology or charges has previously been passed to a calculation method belonging to this class.

Shape: \((N_\mathrm{g},)\).

Reference unit: \(\mathrm{e}\).

calculate_transport_coefficients(start: int = 1, stop: int = None, scale: str = 'log', *, start_self: int = None, stop_self: int = None, scale_self: str = None, enforce_linear: bool = True) None[source]

Fits the mean squared displacements (MSDs) (or analogs) to evaluate the Onsager transport coefficients \(L_{ij}\) and \(L_{ii}^\mathrm{self}\).

Parameters:
startint, default: 1

Starting frame with respect to the interval used in Onsager.run() for fitting the MSDs (or their analogs).

stopint, optional

Ending frame with respect to the interval used in Onsager.run() for fitting the MSDs (or their analogs).

scalestr, {"log", "linear"}

Data scaling for fitting the MSDs (or their analogs) against time.

Valid values:

  • "linear": Linear \(x\)- and \(y\)-axes.

  • "log": Logarithmic \(x\)- and \(y\)-axes.

start_selfint, keyword-only, optional

Starting frame with respect to the interval used in Onsager.run() for fitting the self MSDs. If not provided, start_self shares a value with start.

stop_selfint, keyword-only, optional

Ending frame with respect to the interval used in Onsager.run() for fitting the self MSDs. If not provided, stop_self shares a value with stop.

scale_selfstr, keyword-only, optional

Data scaling for fitting the self MSDs against time. If not provided, scale_self shares a value with scale.

Valid values:

  • "linear": Linear \(x\)- and \(y\)-axes.

  • "log": Logarithmic \(x\)- and \(y\)-axes.

enforce_linearbool, keyword-only, default: True

Enforce linear fits for data on a logarithmic scale by setting the slope to \(1\), or

\[\log(\mathrm{MSD})=\log(t)+b\]
run(start: int = None, stop: int = None, step: int = None, frames: slice | ndarray[int] = None, n_jobs: int = 1, verbose: bool = None, **kwargs) SerialAnalysisBase

Performs the calculation in serial.

Parameters:
startint, optional

Starting frame for analysis.

stopint, optional

Ending frame for analysis.

stepint, optional

Number of frames to skip between each analyzed frame.

framesslice or array-like, optional

Index or logical array of the desired trajectory frames.

verbosebool, optional

Determines whether detailed progress is shown.

**kwargs

Additional keyword arguments to pass to MDAnalysis.lib.log.ProgressBar.

Returns:
selfSerialAnalysisBase

Serial analysis base object.

save(file: str | TextIO, archive: bool = True, compress: bool = True, **kwargs) None

Saves results to a binary or archive file in NumPy format.

Parameters:
filestr or file

Filename or file-like object where the data will be saved. If file is a str, the .npy or .npz extension will be appended automatically if not already present.

archivebool, default: True

Determines whether the results are saved to a single archive file. If True, the data is stored in a .npz file. Otherwise, the data is saved to multiple .npy files.

compressbool, default: True

Determines whether the .npz file is compressed. Has no effect when archive=False.

**kwargs

Additional keyword arguments to pass to numpy.save(), numpy.savez(), or numpy.savez_compressed(), depending on the values of archive and compress.