Anisotropic Network Model¶
This module defines a class and a function for anisotropic network model (ANM) calculations.
-
class
ANM(name='Unknown')[source]¶ Class for Anisotropic Network Model (ANM) analysis of proteins ([PD00], [ARA01]).
See a usage example in Anisotropic Network Model (ANM).
[PD00] Doruker P, Atilgan AR, Bahar I. Dynamics of proteins predicted by molecular dynamics simulations and analytical approaches: Application to a-amylase inhibitor. Proteins 2000 40:512-524. [ARA01] Atilgan AR, Durrell SR, Jernigan RL, Demirel MC, Keskin O, Bahar I. Anisotropy of fluctuation dynamics of proteins with an elastic network model. Biophys. J. 2001 80:505-515. -
addEigenpair(vector, value=None)¶ Add eigen vector and eigen value pair(s) to the instance. If eigen value is omitted, it will be set to 1. Inverse eigenvalues are set as variances.
-
buildHessian(coords, cutoff=15.0, gamma=1.0, **kwargs)¶ Build Hessian matrix for given coordinate set.
Parameters: - coords (
numpy.ndarray) – a coordinate set or an object withgetCoordsmethod - cutoff (float) – cutoff distance (Å) for pairwise interactions, default is 15.0 Å, minimum is 4.0 Å
- gamma (float,
Gamma) – spring constant, default is 1.0 - sparse (bool) – elect to use sparse matrices, default is False. If
Scipy is not found,
ImportErroris raised. - kdtree (bool) – elect to use KDTree for building Hessian matrix, default is False since KDTree method is slower
Instances of
Gammaclasses and custom functions are accepted as gamma argument.When Scipy is available, user can select to use sparse matrices for efficient usage of memory at the cost of computation speed.
- coords (
-
buildMechStiff(coords, n_modes=None, kbt=1.0)¶ Calculate stiffness matrix calculated using
ANMinstance. Method described in [EB08].[EB08] Eyal E., Bahar I. Toward a Molecular Understanding of the Anisotropic Response of Proteins to External Forces: Insights from Elastic Network Models. Biophys J 2008 94:3424-34355. Parameters: - coords (
numpy.ndarray.) – a coordinate set or an object withgetCoordsmethod - n_modes (int or
None, default is 20.) – number of non-zero eigenvalues/vectors to calculate. IfNoneis given, all modes will be calculated (3x number of atoms).
Author: Mustafa Tekpinar & Karolina Mikulska-Ruminska & Cihan Kaya
- coords (
-
calcModes(n_modes=20, zeros=False, turbo=True)¶ Calculate normal modes. This method uses
scipy.linalg.eigh()function to diagonalize the Hessian matrix. When Scipy is not found,numpy.linalg.eigh()is used.Parameters: - n_modes (int or None, default is 20) – number of non-zero eigenvalues/vectors to calculate.
If
Noneis given, all modes will be calculated. - zeros (bool, default is
False) – IfTrue, modes with zero eigenvalues will be kept. - turbo (bool, default is
True) – Use a memory intensive, but faster way to calculate modes.
- n_modes (int or None, default is 20) – number of non-zero eigenvalues/vectors to calculate.
If
-
getArray()¶ Returns a copy of eigenvectors array.
-
getCovariance()¶ Returns covariance matrix. If covariance matrix is not set or yet calculated, it will be calculated using available modes.
-
getCutoff()¶ Returns cutoff distance.
-
getEigvals()¶ Returns eigenvalues. For
PCAandEDAmodels built using coordinate data in Å, unit of eigenvalues is Å2. ForANM,GNM, andRTB, on the other hand, eigenvalues are in arbitrary or relative units but they correlate with stiffness of the motion along associated eigenvector.
-
getEigvecs()¶ Returns a copy of eigenvectors array.
-
getGamma()¶ Returns spring constant (or the gamma function or
Gammainstance).
-
getHessian()¶ Returns a copy of the Hessian matrix.
-
getKirchhoff()¶ Returns a copy of the Kirchhoff matrix.
-
getMechStiffStatistic(rangeK, minAA=0, AA='all')¶ Returns number of effective spring constant with set range of amino acids of protein structure.
AAcan be a list with a range of analysed amino acids as: [first_aa, last_aa, first_aa2, last_aa2], minAA - eliminate amino acids that are within 20aa andrangeKis a list [minK, maxK]
-
getModel()¶ Returns self.
-
getStiffness()¶ Return a copy of Stiffness matrix.
-
getStiffnessRange()¶ Return the range of effective spring constant.
-
getStiffnessRangeSel(value, minAA=20, AA='all')¶ Returns minimum or maximum value of sping constant from mechanical stiffness calculations for residues that are within more than
min_aafrom each other.Valueshould be ‘minK’ or ‘maxK’. It alow to avoid residues near each other.AAis a number of residues from both terminus (N and C) of protein strcuture, it can beallor int value (than first and lastAAresidues will be analysed. WithminAA=0it can be used to search the highest/lowest values of interactions between N-C terminus if protein structure has a shear, zipper or SD1-disconnected mechanical clamp -it is common in FnIII/Ig like domains and determines the maximum unfolding force in AFM or SMD method.
-
getTitle()¶ Returns title of the model.
-
getVariances()¶ Returns variances. For
PCAandEDAmodels built using coordinate data in Å, unit of variance is Å2. ForANM,GNM, andRTB, on the other hand, variance is the inverse of the eigenvalue, so it has arbitrary or relative units.
-
is3d()¶ Returns True if model is 3-dimensional.
-
numAtoms()¶ Returns number of atoms.
-
numDOF()¶ Returns number of degrees of freedom.
-
numModes()¶ Returns number of modes in the instance (not necessarily maximum number of possible modes).
-
setEigens(vectors, values=None)¶ Set eigen vectors and eigen values. If eigen values are omitted, they will be set to 1. Inverse eigenvalues are set as variances.
-
setHessian(hessian)¶ Set Hessian matrix. A symmetric matrix is expected, i.e. not a lower- or upper-triangular matrix.
-
setTitle(title)¶ Set title of the model.
-