Residence

Residence

Author:

Orion Cohen, Tingzheng Hou, Kara Fong

Year:

2021

Copyright:

GNU Public License v3

Understand the dynamic coordination of solvents with the solute.

Residence times for all solvents are automatically calculated from autocovariance of the solvent-solute adjacency matrix.

While residence can be used in isolation, it is meant to be used as an attribute of the Solute class. This makes instantiating it and calculating the solvation data a non-issue.

class solvation_analysis.residence.Residence[source]

Calculate the residence times of solvents.

This class calculates the residence time of each solvent on the solute. The residence time is in units of Solute frames, so if the Solute object has 1000 frames over 1 nanosecond, then each frame will be 1 picosecond. Thus a residence time of 100 would translate to 100 picoseconds.

Two residence time implementations are available. Both calculate the solute-solvent autocorrelation function for each solute-solvent pair, take and take the mean over the solvents of each type, this should yield an exponentially decaying autocorrelation function.

The first implementation fits an exponential curve to the autocorrelation function and extract the time constant, which is inversely proportional to the residence time. This result is saved in the residence_times_fit attribute. Unfortunately, the fit often fails to converge (value is set to np.nan), making this method unreliable.

Instead, the default implementation is to simply find point where the value of the autocorrelation function is 1/e, which is the time constant of an exact exponential. These values are saved in residence_times.

It is recommended that the user visually inspect the autocorrelation function with Residence.plot_autocorrelation_function to ensure an approximately exponential decay. The residence times are only valid insofar as the autocorrelation function resembles an exact exponential, it should decays to zero with a long tail. If the exponential does not decay to zero or its slope does not level off, increasing the simulation time may help. For this technique to be appropriate, the simulation time should exceed the residence time.

A fuller description of the method can be found in Self, Fong, and Persson

Parameters:
  • solvation_data (pandas.DataFrame) – The solvation data frame output by Solute.

  • step (int) – The spacing of frames in solvation_data. This should be equal to solute.step.

Examples

# first define Li, BN, and FEC AtomGroups
>>> solute = Solute(Li, {'BN': BN, 'FEC': FEC, 'PF6': PF6})
>>> residence = Residence.from_solute(solute)
>>> residence.residence_times_cutoff
{'BN': 4.02, 'FEC': 3.79, 'PF6': 1.15}
__init__(solvation_data, step)[source]
static from_solute(solute)[source]

Generate a Residence object from a solute.

Parameters:

solute (Solute)

Return type:

Residence

plot_auto_covariance(res_name)[source]

Plot the autocovariance of a solvent on the solute.

See the discussion in the class docstring for more information.

Parameters:

res_name (str) – the name of a solvent in the solute.

Returns:

  • fig (matplotlib.Figure)

  • ax (matplotlib.Axes)

property auto_covariances

A dictionary where keys are residue names and values are the autocovariance of the that residue on the solute.

property residence_times_cutoff

A dictionary where keys are residue names and values are the residence times of the that residue on the solute, calculated with the 1/e cutoff method.

property residence_times_fit

A dictionary where keys are residue names and values are the residence times of the that residue on the solute, calculated with the exponential fit method.

property fit_parameters

A dictionary where keys are residue names and values are the arameters for the exponential fit to the autocorrelation function.