4.8.1.2. HELANAL — analysis of protein helices

Author:

Lily Wang

Year:

2020

Copyright:

GNU Public License v3

Added in version 2.0.0.

This module contains code to analyse protein helices using the HELANAL algorithm ([Bansal2000] , [Sugeta1967] ).

HELANAL quantifies the geometry of helices in proteins on the basis of their Cα atoms. It can determine local structural features such as the local helical twist and rise, virtual torsion angle, local helix origins and bending angles between successive local helix axes.

[Sugeta1967]

Sugeta, H. and Miyazawa, T. 1967. General method for calculating helical parameters of polymer chains from bond lengths, bond angles and internal rotation angles. Biopolymers 5 673 - 679

[Bansal2000]

Bansal M, Kumar S, Velavan R. 2000. HELANAL - A program to characterise helix geometry in proteins. J Biomol Struct Dyn. 17(5):811-819.

4.8.1.2.1. Example use

You can pass in a single selection:

import MDAnalysis as mda
from MDAnalysis.tests.datafiles import PSF, DCD
from MDAnalysis.analysis import helix_analysis as hel
u = mda.Universe(PSF, DCD)
helanal = hel.HELANAL(u, select='name CA and resnum 161-187')
helanal.run()

All computed properties are available in .results:

print(helanal.results.summary)

Alternatively, you can analyse several helices at once by passing in multiple selection strings:

helanal2 = hel.HELANAL(u, select=('name CA and resnum 100-160',
                                  'name CA and resnum 200-230'))

The helix_analysis() function will carry out helix analysis on atom positions, treating each row of coordinates as an alpha-carbon equivalent:

hel_xyz = hel.helix_analysis(u.atoms.positions, ref_axis=[0, 0, 1])

4.8.1.2.2. Classes

class MDAnalysis.analysis.helix_analysis.HELANAL(universe, select='name CA', ref_axis=[0, 0, 1], verbose=False, flatten_single_helix=True, split_residue_sequences=True)[source]

Perform HELANAL helix analysis on your trajectory.

Parameters:
  • universe (Universe or AtomGroup) – The Universe or AtomGroup to apply the analysis to.

  • select (str or iterable of str, optional) – The selection string to create an atom selection that the HELANAL analysis is applied to. Note that HELANAL is designed to work on the alpha-carbon atoms of protein residues. If you pass in multiple selections, the selections will be analysed separately.

  • ref_axis (array-like of length 3, optional) – The reference axis used to calculate the tilt of the vector of best fit, and the local screw angles.

  • flatten_single_helix (bool, optional) – Whether to flatten results if only one selection is passed.

  • split_residue_sequences (bool, optional) – Whether to split the residue sequence into individual helices. This keyword only applies if a residue gap is present in the AtomGroup generated by a select string. If False, the residues will be analysed as a single helix. If True, each group of consecutive residues will be treated as a separate helix.

  • verbose (bool, optional) – Turn on more logging and debugging.

results.local_twists

The local twist angle from atom i+1 to i+2. Each array has shape (n_frames, n_residues-3)

Type:

array or list of arrays

results.local_nres_per_turn

Number of residues per turn, based on local_twist. Each array has shape (n_frames, n_residues-3)

Type:

array or list of arrays

results.local_axes

The length-wise helix axis of the local window. Each array has shape (n_frames, n_residues-3, 3)

Type:

array or list of arrays

results.local_heights

The rise of each local helix. Each array has shape (n_frames, n_residues-3)

Type:

array or list of arrays

results.local_helix_directions

The unit vector from each local origin to atom i+1. Each array has shape (n_frames, n_residues-2, 3)

Type:

array or list of arrays

results.local_origins

The projected origin for each helix. Each array has shape (n_frames, n_residues-2, 3)

Type:

array or list of arrays

results.local_screw_angles

The local screw angle for each helix. Each array has shape (n_frames, n_residues-2)

Type:

array or list of arrays

results.local_bends

The angles between local helix axes, 3 windows apart. Each array has shape (n_frames, n_residues-6)

Type:

array or list of arrays

results.all_bends

The angles between local helix axes. Each array has shape (n_frames, n_residues-3, n_residues-3)

Type:

array or list of arrays

results.global_axis

The length-wise axis for the overall helix. This points at the first helix window in the helix, so it runs opposite to the direction of the residue numbers. Each array has shape (n_frames, 3)

Type:

array or list of arrays

results.global_tilts

The angle between the global axis and the reference axis. Each array has shape (n_frames,)

Type:

array or list of arrays

results.summary

Summary of stats for each property: the mean, the sample standard deviation, and the mean absolute deviation.

Type:

dict or list of dicts

4.8.1.2.3. Functions

MDAnalysis.analysis.helix_analysis.helix_analysis(positions, ref_axis=[0, 0, 1])[source]

Calculate helix properties from atomic coordinates.

Each property is calculated from a sliding window of 4 atoms, from i to i+3. Any property whose name begins with ‘local’ is a property of a sliding window.

Parameters:
  • positions (numpy.ndarray of shape (N, 3)) – Atomic coordinates.

  • ref_axis (array-like of length 3, optional) – The reference axis used to calculate the tilt of the vector of best fit, and the local screw angles.

Returns:

local_twistsarray, shape (N-3,)

local twist angle from atom i+1 to i+2

local_nres_per_turnarray, shape (N-3,)

number of residues per turn, based on local_twist

local_axesarray, shape (N-3, 3)

the length-wise helix axis of the local window

local_bendsarray, shape (N-6,)

the angles between local helix angles, 3 windows apart

local_heightsarray, shape (N-3,)

the rise of each local helix

local_helix_directionsarray, shape (N-2, 3)

the unit vector from each local origin to atom i+1

local_originsarray, shape (N-2, 3)

the projected origin for each helix

all_bendsarray, shape (N-3, N-3)

angles between each local axis

global_axisarray, shape (3,)

vector of best fit through origins, pointing at the first origin.

local_screw_anglesarray, shape (N-2,)

cylindrical azimuth angle to plane of global_axis and ref_axis

Return type:

dict with the following keys

MDAnalysis.analysis.helix_analysis.vector_of_best_fit(coordinates)[source]

Fit vector through the centered coordinates, pointing to the first coordinate (i.e. upside-down).

Parameters:

coordinates (numpy.ndarray of shape (N, 3))

Returns:

Vector of best fit.

Return type:

numpy.ndarray of shape (3,)

MDAnalysis.analysis.helix_analysis.local_screw_angles(global_axis, ref_axis, helix_directions)[source]

Cylindrical azimuth angles between the local direction vectors, as projected onto the cross-section of the helix, from (-pi, pi]. The origin (angle=0) is set to the plane of global_axis and ref_axis.

Parameters:
  • global_axis (numpy.ndarray of shape (3,)) – Vector of best fit. Screw angles are calculated perpendicular to this axis.

  • ref_axis (numpy.ndarray of shape (3,)) – Reference length-wise axis. One of the reference vectors is orthogonal to this axis.

  • helix_directions (numpy.ndarray of shape (N, 3)) – array of vectors representing the local direction of each helix window.

Returns:

Array of screw angles.

Return type:

numpy.ndarray of shape (N,)