Explicit Membrane Gaussian Network Model

This module defines a class and a function for explicit membrane GNM calculations.

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

Class for explicit GNM (exGNM) method ([FT00]). Optional arguments build a membrane lattice permit analysis of membrane

effect on elastic network models in exGNM method described in [TL12].
[TL12]Lezon TR, Bahar I, Constraints Imposed by the Membrane Selectively Guide the Alternating Access Dynamics of the Glutamate Transporter GltPh
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=15.0, gamma=1.0, **kwargs)[source]

Build Kirchhoff matrix for given coordinate set. kwargs are passed to :method:`.buildMembrane`.

Parameters:
  • coords (numpy.ndarray) – a coordinate set or an object with getCoords method
  • cutoff (float) – cutoff distance (Å) for pairwise interactions, default is 15.0 Å
  • gamma (float) – spring constant, default is 1.0
buildMembrane(coords, **kwargs)[source]

Build Kirchhoff matrix for given coordinate set.

Parameters:
  • coords (numpy.ndarray) – a coordinate set or an object with getCoords method
  • membrane_high (float) – the maximum z coordinate of the membrane. Default is 13.0
  • membrane_low (float) – the minimum z coordinate of the membrane. Default is -13.0
  • R (float) – radius of all membrane in x-y direction. Default is 80
  • Ri (float) – inner radius of the membrane in x-y direction if it needs to be hollow. Default is 0, which is not hollow
  • r (float) – radius of each membrane node. Default is .1*
  • lat (str) – lattice type which could be FCC (face-centered-cubic, default), SC (simple cubic), SH (simple hexagonal)
  • exr (float) – exclusive radius of each protein node. Default is 5.0
  • hull (bool) – whether use convex hull to determine the protein’s interior. Turn it off if protein is multimer. Default is True
  • center (bool) – whether transform the structure to the origin (only x- and y-axis). Default is True
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 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.

getCombined()[source]

Returns a copy of the combined atoms or coordinates.

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.

getMembrane()[source]

Returns a copy of the membrane coordinates.

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

Set Kirchhoff matrix.

setTitle(title)

Set title of the model.