Analysis Functions

This module defines functions for calculating physical properties from normal modes.

calcCollectivity(mode, masses=None)[source]

Returns collectivity of the mode. This function implements collectivity as defined in equation 5 of [BR95]. If masses are provided, they will be incorporated in the calculation. Otherwise, atoms are assumed to have uniform masses.

[BR95]Bruschweiler R. Collective protein dynamics and nuclear spin relaxation. J Chem Phys 1995 102:3396-3403.
Parameters:
calcCovariance(modes)[source]

Returns covariance matrix calculated for given modes.

calcCrossCorr(modes, n_cpu=1)[source]

Returns cross-correlations matrix. For a 3-d model, cross-correlations matrix is an NxN matrix, where N is the number of atoms. Each element of this matrix is the trace of the submatrix corresponding to a pair of atoms. Covariance matrix may be calculated using all modes or a subset of modes of an NMA instance. For large systems, calculation of cross-correlations matrix may be time consuming. Optionally, multiple processors may be employed to perform calculations by passing n_cpu=2 or more.

calcFractVariance(mode)[source]

Returns fraction of variance explained by the mode. Fraction of variance is the ratio of the variance along a mode to the trace of the covariance matrix of the model.

calcSqFlucts(modes)[source]

Returns sum of square-fluctuations for given set of normal modes. Square fluctuations for a single mode is obtained by multiplying the square of the mode array with the variance (Mode.getVariance()) along the mode. For PCA and EDA models built using coordinate data in Å, unit of square-fluctuations is Å2, for ANM and GNM, on the other hand, it is arbitrary or relative units.

calcTempFactors(modes, atoms)[source]

Returns temperature (β) factors calculated using modes from a ANM or GNM instance scaled according to the experimental β-factors from atoms.

calcProjection(ensemble, modes, rmsd=True)[source]

Returns projection of conformational deviations onto given modes. ensemble coordinates are used to calculate the deviations that are projected onto modes. For K conformations and M modes, a (K,M) matrix is returned.

Parameters:

By default root-mean-square deviation (RMSD) along the normal mode is calculated. To calculate the projection pass rmsd=True. Vector instances are accepted as ensemble argument to allow for projecting a deformation vector onto normal modes.

calcCrossProjection(ensemble, mode1, mode2, scale=None, **kwargs)[source]

Returns projection of conformational deviations onto modes from different models.

Parameters:
  • ensemble (Ensemble) – ensemble for which deviations will be projected
  • mode1 (Mode, Vector) – normal mode to project conformations onto
  • mode2 (Mode, Vector) – normal mode to project conformations onto
  • scale – scale width of the projection onto mode1 (x) or mode2(y), an optimized scaling factor (scalar) will be calculated by default or a value of scalar can be passed.
calcPerturbResponse(model, atoms=None, repeats=100, **kwargs)[source]

Returns a matrix of profiles from scanning of the response of the structure to random perturbations at specific atom (or node) positions. The function implements the perturbation response scanning (PRS) method described in [CA09]. Rows of the matrix are the average magnitude of the responses obtained by perturbing the atom/node position at that row index, i.e. prs_profile[i,j] will give the response of residue/node j to perturbations in residue/node i. PRS is performed using the covariance matrix from model, e.g. ANM instance. Each residue/node is perturbed repeats times with a random unit force vector. When atoms instance is given, PRS profile for residues will be added as an attribute which then can be retrieved as atoms.getData('prs_profile'). model and atoms must have the same number of atoms. atoms must be an AtomGroup instance. write_output is a Bool variable to write normalized asymmetric PRS matrix to a file.

[CA09]Atilgan C, Atilgan AR, Perturbation-Response Scanning Reveals Ligand Entry-Exit Mechanisms of Ferric Binding Protein. PLoS Comput Biol 2009 5(10):e1000544.

The PRS matrix can be calculated and saved as follows:

prs_matrix = calcPerturbationResponse(p38_anm, saveMatrix=True)

The PRS matrix can also be save later as follows:

writeArray('prs_matrix.txt', prs_matrix, format='%8.6f', delimiter='      ')

You can also control which operation is used for getting a single matrix from the repeated force application and whether to normalise the matrix at the end. If you do choose to normalise the matrix, you can still save the original matrix before normalisation as well.

Parameters:
  • operation (str or list) – which operation to perform to get a single response matrix:: the mean, variance, max or min of the set of repeats. Another operation is to select elements from the matrix showing biggest difference from the square sum of the covariance matrix. The Default is the mean. To obtain all response matrices, set operation=None without quotes. You can also ask for ‘all’ operations or provide a list containing any set of them.
  • noForce (bool) – whether to use the covariance matrix directly rather than applying forces. This appears to be equivalent when scanning for response magnitudes and will be much quicker. Default is True.
  • normMatrix (bool) – whether to normalise the single response matrix by dividing each row by its diagonal, Default is False, we recommend true
  • saveMatrix (bool) – whether to save the last matrix generated to a text file. Default is False
  • saveOrig (bool) – whether to save the original matrix despite normalisation. This is the same as saveMatrix when not normalizing. Default is False
  • baseSaveName (str) – The central part of the file name for saved matrices, which you can set. This is surrounded by underscores. The beginning says orig or norm and the end says which operation was used. Default is ‘response_matrix’.
  • acceptDirection (str) – select reference direction for forces to be accepted. Can be ‘in’ (towards center of atoms), ‘out’ (away from center), or ‘all’. Default is ‘all’; Using other directions requires atoms.
parsePerturbResponseMatrix(prs_matrix_file='prs_matrix.txt', normMatrix=True)[source]

Parses a perturbation response matrix from a file into a numpy ndarray.

Parameters:
  • prs_matrix_file (str) – name of the file containing a PRS matrix, default is ‘prs_matrix.txt’ as is used in the example under calcPerturbResponse.
  • normMatrix – whether to normalise the PRS matrix after parsing it. Default is True. If you have already normalised this before saving, or you do not want it normalised then set this to False.
calcPerturbResponseProfiles(prs_matrix)[source]

Calculate the effectiveness and sensitivity profiles, which are the averages over the rows and columns of the PRS matrix.

Parameters:prs_matrix (ndarray) – a perturbation response matrix
writePerturbResponsePDB(prs_matrix, pdbIn, **kwargs)[source]

Write the average response to perturbation of a particular residue (a row of a perturbation response matrix) or the average effect of perturbation of a particular residue (a column of a normalized perturbation response matrix) into the b-factor field of a PDB file for visualisation in a molecular graphics program. If no chain is given this will be done for that residue in all chains.

If no residue number is given then the effectiveness and sensitivity profiles will be written out instead. These two profiles are also returned as arrays for further analysis if they aren’t already provided.

Parameters:
  • prs_matrix (ndarray) – a perturbation response matrix
  • pdbIn (str) – file name for the input PDB file where you would like the PRS data mapped
  • pdbOut (list) – a list of file names (enclosed in square brackets) for the output PDB file, default is to append the chain and residue info (name and number) onto the pdbIn stem. The input for pdbOut can also be used as a stem if you enter a single string enclosed in quotes. If no residue number is supplied, chain is ignored and the default is to append ‘_effectiveness’ and ‘_sensitivity’ onto the stem.
  • chain (str) – chain identifier for the residue of interest, default is all chains If you want to analyse residues in a subset of chains, concatentate them together e.g. ‘AC’
  • resnum (int) – residue number for the residue of interest
  • direction (str) – the direction you want to use to read data out of the PRS matrix for plotting: the options are ‘row’ or ‘column’. Default is ‘row’. A row gives the effect on each residue of peturbing the specified residue. A column gives the response of the specified residue to perturbing each residue. If no residue number is provided then this option will be ignored
  • returnData – whether to return effectiveness and sensitivity for analysis default is False
calcSpecDimension(mode)[source]
Parameters:mode (Mode or Vector) – mode or vector
calcPairDeformationDist(model, coords, ind1, ind2, kbt=1.0)[source]

Returns distribution of the deformations in the distance contributed by each mode for selected pair of residues ind1 ind2 using model from a ANM. Method described in [EB08] equation (10) and figure (2).

[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:model – this is an 3-dimensional NMA instance from a :class:`.ANM

calculations. :type model: ANM :arg coords: a coordinate set or an object with getCoords method.

Recommended: coords = parsePDB(‘pdbfile’).select(‘protein and name CA’).
Parameters:
  • ind1 (int) – first residue number.
  • ind2 (int) – secound residue number.