Gaussian Network Model

This module defines a class and a function for Gaussian network model (GNM) calculations.

class GNM(name='Unknown')[source]

A class for Gaussian Network Model (GNM) analysis of proteins ([IB97], [TH97]).

See example Gaussian Network Model (GNM).

[IB97]Bahar I, Atilgan AR, Erman B. Direct evaluation of thermal fluctuations in protein using a single parameter harmonic potential. Folding & Design 1997 2:173-181.
[TH97]Haliloglu T, Bahar I, Erman B. Gaussian dynamics of folded proteins. Phys. Rev. Lett. 1997 79:3090-3093.
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.

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

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 Kirchhoff 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.
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 Å2. 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).

getKirchhoff()

Returns a copy of the Kirchhoff matrix.

getModel()

Returns self.

getNormDistFluct(coords)[source]

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 Å2. 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).

setKirchhoff(kirchhoff)[source]

Set Kirchhoff matrix.

setTitle(title)

Set title of the model.

calcGNM(pdb, selstr='calpha', cutoff=10, gamma=1.0, n_modes=20, zeros=False)[source]

Returns a GNM 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 Atomic instance.