pylipid.api.LipidInteraction.compute_site_koff

LipidInteraction.compute_site_koff(binding_site_id=None, nbootstrap=10, initial_guess=[1.0, 1.0, 1.0, 1.0], save_dir=None, plot_data=True, fig_close=True, fig_format='pdf', num_cpus=None)[source]

Calculate interactions koff and residence time for binding sites.

The koff is calculated from a survival time correlation function which describes the relaxation of the bound lipids 1. Often the interactions between lipid and protein surface are be divided into prolonged interactions and quick diffusive contacts. Thus PyLipID fits the normalised survival function to a bi-exponential curve which describes the long and short decay periods.

The survival time correlation function σ(t) is calculated as follow

\[\sigma(t) = \frac{1}{N_{j}} \frac{1}{T-t} \sum_{j=1}^{N_{j}} \sum_{v=0}^{T-t}\tilde{n}_{j}(v, v+t)\]

where T is the length of the simulation trajectory, \(N_{j}\) is the total number of lipid contacts and \(\sum_{v=0}^{T-t} \tilde{n}_{j}(v, v+t)\) is a binary function that takes the value 1 if the contact of lipid j lasts from time ν to time v+t and 0 otherwise. The values of \(\sigma(t)\) are calculated for every value of t from 0 to T ns, for each time step of the trajectories, and normalized by dividing by \(\sigma(t)\), so that the survival time-correlation function has value 1 at t = 0.

The normalized survival function is then fitted to a biexponential to model the long and short decays of lipid relaxation:

\[\sigma(t) \sim A e^{-k_{1} t}+B e^{-k_{2} t}\left(k_{1} \leq k_{2}\right)\]

PyLipID takes \(k_{1}\) as the the dissociation \(k_{off}\), and calculates the residence time as \(\tau=1 / k_{o f f}\). PyLipID raises a warning for the impact on the accuracy of \(k_{off}\) calculation if trajectories are of different lengths when multiple trajectories are provided. PyLipID measures the \(r^{2}\) of the biexponential fitting to the survival function to show the quality of the \(k_{off}\) estimation. In addition, PyLipID bootstraps the contact durations and measures the \(k_{off}\) of the bootstrapped data, to report how well lipid contacts are sampled from simulations. The lipid contact sampling, the curve-fitting and the bootstrap results can be conveniently checked via the \(k_{off}\) plot.

The durations of lipid contact with binding sites are calculated using compute_site_duration(). See its page for the definition of lipid contact with binding site.

The calculation of koff for binding sites can be time-consuming, thus PyLipID uses python multiprocessing to parallize the calculation. The number of CPUs used for multiprocessing can be specificed, otherwise all the available CPUs will be used by default.

Parameters
  • binding_site_id (int or list of int, default=None) – The binding site ID used in PyLipID. This ID is the index in the binding site node list that is calculated by the method compute_binding_nodes. The ID of the N-th binding site is (N-1). If None, the contact duration of all binding sites are calculated.

  • nbootstrap (int, default=10) – Number of bootstrap on the interaction durations. For each bootstrap, samples of the size of the original dataset are drawn from the collected durations with replacement. \(k_{koff}\) and \(r^{2}\) are calculated for each bootstrap.

  • initial_guess (array_like, default=None) – The initial guess for the curve-fitting of the biexponential curve. Used by scipy.optimize.curve_fit.

  • save_dir (str, default=None) – The the directory for saving the koff figures of residues if plot_data is True. By default, the koff figures are saved in the directory of Binding_Sites_koffs_{lipid} under the root directory defined when LipidInteraction was initiated.

  • plot_data (bool, default=True) – If True, plot the koff figures fir residues.

  • fig_close (bool, default=True) – Use matplotlib.pyplot.close() to close the koff figures. Can save memory if many figures are open and plotted.

  • fig_format (str, default="pdf") – The format of koff figures. Support formats that are supported by matplotlib.pyplot.savefig().

  • num_cpus (int or None, default=None) – Number of CPUs used for multiprocessing. If None, all the available CPUs will be used.

Returns

  • koff (scalar or list of scalar) – The calculated koffs for selected binding sites.

  • restime (scalar or list of scalar) – The calculated residence times for selected binding sites.

See also

pylipid.api.LipidInteraction.collect_residue_contacts

Create the lipid index.

pylipid.api.LipidInteraction.compute_residue_koff

Calculate koffs and residence times for residues.

pylipid.func.cal_koff

Calculate residence time and koff.

pylipid.func.cal_survival_func

Compute the normalised survival function.

References

1

García, Angel E.Stiller, Lewis. Computation of the mean residence time of water in the hydration shells of biomolecules. 1993. Journal of Computational Chemistry.