Feature of Molecular System — molann.feature
- Author:
Wei Zhang
- Year:
2022
- Copyright:
GNU Public License v3
This module implements a class that defines a feature of a molecular system
(molann.feature.Feature
), and a class that constructs a list of
features from a feature file (read by molann.feature.FeatureFileReader
).
Classes
- class molann.feature.Feature(name, feature_type, atom_group)[source]
Feature of system
- Parameters:
name (str) – feature’s name
feature_type (str) – currently supported values: ‘angle’, ‘bond’, ‘dihedral’, and ‘position’
atom_group (
MDAnalysis.core.groups.AtomGroup
) – atom group used to define a feature
- Raises:
AssertionError – if the number of atoms in atom_group does not match feature_type.
Note
For feature_type = ‘angle’, ‘bond’, and ‘dihedral’, atom_group must contain 3 atoms, 2 atoms, and 4 atoms, respectively.
Examples:
# package MDAnalysis is required import MDAnalysis as mda from molann.feature import Feature # pdb file of the system pdb_filename = '/path/to/system.pdb' ref = mda.Universe(pdb_filename) # feature that is the angle formed by the atoms whose (1-based) indices are 1, 2 and 3. f1 = Feature('some name', 'angle', ref.select_atoms('bynum 1 2 3')) print (f1.get_feature_info()) # feature that is the bond distance between the two atoms whose (1-based) indices are 1 and 2. f2 = Feature('some name', 'bond', ref.select_atoms('bynum 1 2')) # dihedral angle formed by atoms whose (1-based) indices are 1, 2, 3, and 4. f3 = Feature('some name', 'dihedral', ref.select_atoms('bynum 1 2 3 4')) # coordinates of atoms whose (1-based) indices are 3 and 5. f4 = Feature('some name', 'position', ref.select_atoms('bynum 3 5'))
Note
MDAnalysis.core.universe.Universe.select_atoms()
returns an atom group that may not preserve the order of atoms. To construct a feature that defines the dihedral angle formed by atoms whose indices are 2,1,3, and 4, we can define the atom group by concatenation:ag = ref.select_atoms('bynum 2') + ref.select_atoms('bynum 1') + ref.select_atoms('bynum 3 and 4') f = Feature('some name', 'dihedral', ag)
- name
string, name of the feature
- type_name
string, whose value is either ‘angle’, ‘bond’, ‘dihedral’, or ‘position’
- type_id
integer, 0 if type_name=’angle’, 1 if type_name=’bond’, 2 if type_name=’dihedral’, 3 if type_name=’position’
- atom_group
MDAnalysis.core.groups.AtomGroup
used to define the feature
- get_atom_indices()[source]
- Returns:
numpy array of int, (1-based) indices of atoms in the atom group.
- get_feature_info()[source]
- Returns:
pandas.DataFrame
, which contains feature’s information
- class molann.feature.FeatureFileReader(feature_file, section_name, universe)[source]
Read features from file
- Parameters:
feature_file (str) – name of the feature file
section_name (str) – name of the section in the file from which to read features
universe (
MDAnalysis.core.universe.Universe
) – universe that defines the system
Note
A feature file is a normal text file.
Each section describes a list of features. A section begins with a line of content ‘[section_name]’, and ends with a line of content ‘[End]’ (case sensitive). The function
read()
constructs a list ofFeature
from the section with the corresponding section_name.Lines that describe a feature contain several fields, seperated by comma. The first and the second fields specify the name and the type_name of the feature, respectively. The remaining fields specify the strings for selection to define an atom group (by concatenation). See
MDAnalysis.core.selection
.Lines starting with ‘#’ are comment lines and they are not processed.
Examples
Below is an example of feature file, named as feature.txt.
# This is a comment line. # Lines that describe a feature contain several fields, seperated by comma. # The first and the second fields specify the name and the # type_name of the feature, respectively. The remaining fields specify # the strings for selection to define an atom group by concatenation. # Note: to keep the order of atoms, use 'bynum 5, bynum 2', instead of 'bynum 5 2' [Preprocessing] #position, type C or type O or type N p1, position, resid 2 [End] [Histogram] d1, dihedral, bynum 5, bynum 7, bynum 9, bynum 15 d2, dihedral, bynum 7, bynum 9, bynum 15, bynum 17 b1, bond, bynum 2 5 b2, bond, bynum 5 6 a1, angle, bynum 20, bynum 19, bynum 21 a2, angle, bynum 16, bynum 15, bynum 17 [End] [Output] d1, dihedral, bynum 5 7 9 15 d2, dihedral, bynum 7 9 15 17 [End]
The following code constructs a list of features from the section ‘Histogram’, and a list of features from section ‘Preprocessing’.
# package MDAnalysis is required import MDAnalysis as mda from molann.feature import FeatureFileReader # pdb file of the system pdb_filename = '/path/to/system.pdb' ref = mda.Universe(pdb_filename) # read features from the section 'Histogram' feature_reader = FeatureFileReader('feature.txt', 'Histogram', ref) feature_list_1 = feature_reader.read() # read features from the section 'Preprocessing' feature_reader = FeatureFileReader('feature.txt', 'Preprocessing', ref) feature_list_2 = feature_reader.read()
- get_feature_info()[source]
- Returns:
pandas.DataFrame
, which contains information of all features (each row describes a feature)