Anisotropic Network Model

This module defines a class and a function for anisotropic network model (ANM) calculations.

class prody.dynamics.anm.ANM(name='Unknown')[source]

Class for Anisotropic Network Model (ANM) analysis of proteins ([PD00], [ARA01]).

See a usage example in 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.

__getitem__(index)

A list or tuple of integers can be used for indexing.

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, Atomic, Ensemble, Trajectory) – a coordinate set or an object with getCoords method

  • 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, ImportError is raised.

  • kdtree (bool) – elect to use KDTree for building Hessian matrix, default is False since KDTree method is slower

Instances of Gamma classes 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.

Any atoms or points can be used for building a Hessian matrix, including calphas, phosphorus and carbon atoms from nucleic acids, all atoms, or pseudoatoms fitted to density maps with algorithms such as TRN.

The cutoff distance may need to be adjusted depending on the coarse graining level of the atoms or points used.

buildKirchhoff(coords, cutoff=10.0, gamma=1.0, **kwargs)

Build Kirchhoff matrix for given coordinate set.

Parameters:
  • coords (numpy.ndarray or Atomic) – a coordinate set or an object with getCoords method

  • cutoff (float) – cutoff distance (Å) for pairwise interactions default is 10.0 Å, , minimum is 4.0 Å

  • gamma (float) – spring constant, default is 1.0

  • sparse (bool) – elect to use sparse matrices, default is False. If Scipy is not found, ImportError is raised.

  • kdtree (bool) – elect to use KDTree for building Kirchhoff matrix faster, default is True

Instances of Gamma classes 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.

calcModes(n_modes=20, zeros=False, turbo=True, **kwargs)

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 None or 'all' is given, all modes will be calculated.

  • zeros (bool, default is True) – If True, modes with zero eigenvalues will be kept.

  • turbo (bool, default is True) – Use a memory intensive, but faster way to calculate modes.

  • nproc (int) – number of processors for thread pool limit, default is 0, meaning don’t impose limit

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 PCA and EDA models built using coordinate data in Å, unit of eigenvalues is |A2|. For ANM, GNM, and RTB, 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 Gamma instance).

getHessian()

Returns a copy of the Hessian matrix.

getKirchhoff()

Returns a copy of the Kirchhoff matrix.

getModel()

Returns self.

getNormDistFluct(coords)

Normalized distance fluctuation

getTitle()

Returns title of the model.

getVariances()

Returns variances. For PCA and EDA models built using coordinate data in Å, unit of variance is |A2|. For ANM, GNM, and RTB, 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.

numEntries()

Returns number of entries in one eigenvector.

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.

setKirchhoff(kirchhoff)

Set Kirchhoff matrix.

setTitle(title)

Set title of the model.

class prody.dynamics.anm.MaskedANM(name='Unknown', mask=False, masked=True)[source]
__getitem__(index)

A list or tuple of integers can be used for indexing.

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, Atomic, Ensemble, Trajectory) – a coordinate set or an object with getCoords method

  • 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, ImportError is raised.

  • kdtree (bool) – elect to use KDTree for building Hessian matrix, default is False since KDTree method is slower

Instances of Gamma classes 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.

Any atoms or points can be used for building a Hessian matrix, including calphas, phosphorus and carbon atoms from nucleic acids, all atoms, or pseudoatoms fitted to density maps with algorithms such as TRN.

The cutoff distance may need to be adjusted depending on the coarse graining level of the atoms or points used.

buildKirchhoff(coords, cutoff=10.0, gamma=1.0, **kwargs)

Build Kirchhoff matrix for given coordinate set.

Parameters:
  • coords (numpy.ndarray or Atomic) – a coordinate set or an object with getCoords method

  • cutoff (float) – cutoff distance (Å) for pairwise interactions default is 10.0 Å, , minimum is 4.0 Å

  • gamma (float) – spring constant, default is 1.0

  • sparse (bool) – elect to use sparse matrices, default is False. If Scipy is not found, ImportError is raised.

  • kdtree (bool) – elect to use KDTree for building Kirchhoff matrix faster, default is True

Instances of Gamma classes 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.

calcModes(n_modes=20, zeros=False, turbo=True)[source]

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 None or 'all' is given, all modes will be calculated.

  • zeros (bool, default is True) – If True, modes with zero eigenvalues will be kept.

  • turbo (bool, default is True) – Use a memory intensive, but faster way to calculate modes.

  • nproc (int) – number of processors for thread pool limit, default is 0, meaning don’t impose limit

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 PCA and EDA models built using coordinate data in Å, unit of eigenvalues is |A2|. For ANM, GNM, and RTB, 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 Gamma instance).

getHessian()

Returns a copy of the Hessian matrix.

getKirchhoff()

Returns a copy of the Kirchhoff matrix.

getModel()

Returns self.

getNormDistFluct(coords)

Normalized distance fluctuation

getTitle()

Returns title of the model.

getVariances()

Returns variances. For PCA and EDA models built using coordinate data in Å, unit of variance is |A2|. For ANM, GNM, and RTB, 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.

numEntries()

Returns number of entries in one eigenvector.

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)[source]

Set Hessian matrix. A symmetric matrix is expected, i.e. not a lower- or upper-triangular matrix.

setKirchhoff(kirchhoff)

Set Kirchhoff matrix.

setNumAtoms(n_atoms)

Fixes the tail of the model. If n_atoms is greater than the original size (number of nodes), then extra hidden nodes will be added to the model, and if not, the model will be trimmed so that the total number of nodes, including the hidden ones, will be equal to the n_atoms. Note that if masked is True, the n_atoms should be the number of visible nodes.

Parameters:

n_atoms (int) – number of atoms of the model

setTitle(title)

Set title of the model.

prody.dynamics.anm.calcANM(pdb, selstr='calpha', cutoff=15.0, gamma=1.0, n_modes=20, zeros=False, title=None)[source]

Returns an ANM instance and atoms used for the calculations. By default only alpha carbons are considered, but selection string helps selecting a subset of it. pdb can be a PDB code, Atomic instance, or a Hessian matrix (ndarray).