7.4. Core Topology Objects — MDAnalysis.core.topologyobjects
¶
The building blocks for MDAnalysis’ description of topology
-
class
MDAnalysis.core.topologyobjects.
Angle
(atoms, is_guessed=False)[source]¶ An angle between three
Atom
instances. Atom 2 is the apex of the angleNew in version 0.8.
Changed in version 0.9.0: Now a subclass of
TopologyObject
; now uses__slots__
and stores atoms inatoms
attribute-
angle
()[source]¶ Returns the angle in degrees of this Angle.
Angle between atoms 0 and 2 with apex at 1:
2 / / 1------0
Note
The numerical precision is typically not better than 4 decimals (and is only tested to 3 decimals).
New in version 0.9.0.
-
value
()¶ Returns the angle in degrees of this Angle.
Angle between atoms 0 and 2 with apex at 1:
2 / / 1------0
Note
The numerical precision is typically not better than 4 decimals (and is only tested to 3 decimals).
New in version 0.9.0.
-
-
class
MDAnalysis.core.topologyobjects.
Bond
(atoms, order=None, is_guessed=False)[source]¶ A bond between two
Atom
instances.Two
Bond
instances can be compared with the==
and!=
operators. A bond is equal to another if the same atom numbers are connected and they have the same bond order. The ordering of the two atom numbers is ignored as is the fact that a bond was guessed.The presence of a particular atom can also be queried:
>>> Atom in Bond
will return either
True
orFalse
.Changed in version 0.9.0: Now a subclass of
TopologyObject
. Changed class to use__slots__
and stores atoms inatoms
attribute.-
value
(pbc=False)¶ Length of the bond.
Changed in version 0.11.0: Added pbc keyword
-
-
class
MDAnalysis.core.topologyobjects.
Dihedral
(atoms, is_guessed=False)[source]¶ Dihedral (dihedral angle) between four
Atom
instances.The dihedral is defined as the angle between the planes formed by Atoms (1, 2, 3) and (2, 3, 4).
New in version 0.8.
Changed in version 0.9.0: Now a subclass of
TopologyObject
; now uses__slots__
and stores atoms inatoms
attribute.Changed in version 0.11.0: Renamed to Dihedral (was Torsion)
-
dihedral
()[source]¶ Calculate the dihedral angle in degrees.
Dihedral angle around axis connecting atoms 1 and 2 (i.e. the angle between the planes spanned by atoms (0,1,2) and (1,2,3)):
3 | 1-----2 / 0
Note
The numerical precision is typically not better than 4 decimals (and is only tested to 3 decimals).
New in version 0.9.0.
-
value
()¶ Calculate the dihedral angle in degrees.
Dihedral angle around axis connecting atoms 1 and 2 (i.e. the angle between the planes spanned by atoms (0,1,2) and (1,2,3)):
3 | 1-----2 / 0
Note
The numerical precision is typically not better than 4 decimals (and is only tested to 3 decimals).
New in version 0.9.0.
-
-
class
MDAnalysis.core.topologyobjects.
ImproperDihedral
(atoms, is_guessed=False)[source]¶ Improper Dihedral (improper dihedral angle) between four
Atom
instances.MDAnalysis treats the improper dihedral angle as the angle between the planes formed by Atoms (1, 2, 3) and (2, 3, 4).
Warning
Definitions of Atom ordering in improper dihedrals can change. Check the definitions here against your software.
New in version 0.9.0.
Changed in version 0.11.0: Renamed to ImproperDihedral (was Improper_Torsion)
-
class
MDAnalysis.core.topologyobjects.
TopologyDict
(members)[source]¶ A customised dictionary designed for sorting the bonds, angles and dihedrals present in a group of atoms.
Usage:
topologydict = TopologyDict(members)
Parameters: *members* – A list of TopologyObject
instancesReturns: topologydict A specialised dictionary of the topology instances passed to it TopologyDicts are also built lazily from a
TopologyGroup.topDict
attribute.The
TopologyDict
collects all the selected topology type from the atoms and categorises them according to the types of the atoms within. ATopologyGroup
containing all of a given bond type can be made by querying with the appropriate key. The keys to theTopologyDict
are a tuple of the atom types that the bond represents and can be viewed using thekeys()
method.For example, from a system containing pure ethanol
>>> td = u.bonds.topDict >>> td.keys() [('C', 'C'), ('C', 'H'), ('O', 'H'), ('C', 'O')] >>> td['C', 'O'] < TopologyGroup containing 912 bonds >
Note
The key for a bond is taken from the type attribute of the atoms.
Getting and setting types of bonds is done smartly, so a C-C-H angle is considered identical to a H-C-C angle.
Duplicate entries are automatically removed upon creation and combination of different Dicts. This means a bond between atoms 1 and 2 will only ever appear once in a dict despite both atoms 1 and 2 having the bond in their
bond
attribute.Two
TopologyDict
instances can be combined using addition and it will not create any duplicate bonds in the process.New in version 0.8.
Changed in version 0.9.0: Changed initialisation to use a list of
TopologyObject
instances instead of list of atoms; now used from withinTopologyGroup
instead of accessed fromAtomGroup
.-
__contains__
(other)[source]¶ Returns boolean on whether a given type exists within this dictionary
For topology groups the key (1,2,3) is considered the same as (3,2,1)
-
-
class
MDAnalysis.core.topologyobjects.
TopologyGroup
(bondlist)[source]¶ A container for a groups of bonds.
All bonds of a certain types can be retrieved from within the
TopologyGroup
by querying with a tuple of types:tg2 = tg.select_bonds([key])
Where key describes the desired bond as a tuple of the involved
Atom
types, as defined by the .type Atom attribute). A list of available keys can be displayed using thetypes()
method.Alternatively, all the bonds which are in a given
AtomGroup
can be extracted usingatomgroup_intersection()
:tg2 = tg.atomgroup_intersection(ag)
This allows the keyword strict to be given, which forces all members of all bonds to be inside the AtomGroup passed to it.
Finally, a TopologyGroup can be sliced similarly to AtomGroups:
tg2 = tg[5:10]
The
bonds()
,angles()
anddihedrals()
methods offer a “shortcut” to the Cython distance calculation functions inMDAnalysis.lib.distances
.TopologyGroups can be combined with TopologyGroups of the same bond type (ie can combine two angle containing TopologyGroups).
New in version 0.8.
Changed in version 0.9.0: Overhauled completely: (1) Added internal
TopologyDict
accessible by thetopDict
attribute. (2)selectBonds()
allows thetopDict
to be queried with tuple of types. (3) Addedatomgroup_intersection()
to allow bonds which are in a givenAtomGroup
to be retrieved.Changed in version 0.10.0: Added
from_indices()
constructor, allowing class to be created from indices. Can now create empty Group. Renameddump_contents()
toto_indices()
Changed in version 0.11.0: Added values method to return the size of each object in this group Deprecated selectBonds method in favour of select_bonds
-
__add__
(other)[source]¶ Combine two TopologyGroups together.
Can combined two TopologyGroup of the same type, or add a single TopologyObject to a TopologyGroup.
-
__getitem__
(item)[source]¶ Returns a particular bond as single object or a subset of this TopologyGroup as another TopologyGroup
Changed in version 0.10.0: Allows indexing via boolean numpy array
-
angles
(result=None, pbc=False)[source]¶ Calculates the angle in radians formed between a bond between atoms 1 and 2 and a bond between atoms 2 & 3
Keywords: - result
allows a predefined results array to be used, note that this will be overwritten
- pbc
apply periodic boundary conditions when calculating angles [
False
] this is important when connecting vectors between atoms might require minimum image convention
Uses cython implementation
Changed in version 0.9.0: Added pbc option (default
False
)
-
atomgroup_intersection
(ag, **kwargs)[source]¶ Retrieve all bonds from within this TopologyGroup that are within the AtomGroup which is passed.
Keywords: - strict
Only retrieve bonds which are completely contained within the AtomGroup. [
False
]
New in version 0.9.0.
-
bonds
(pbc=False, result=None)[source]¶ Calculates the distance between all bonds in this TopologyGroup
Keywords: - pbc
apply periodic boundary conditions when calculating distance [False]
- result
allows a predefined results array to be used, note that this will be overwritten
Uses cython implementation
-
dihedrals
(result=None, pbc=False)[source]¶ - Calculate the dihedralal angle in radians for this topology
group.
Defined as the angle between a plane formed by atoms 1, 2 and 3 and a plane formed by atoms 2, 3 and 4.
Keywords: - result
allows a predefined results array to be used, note that this will be overwritten
pbc
- v apply periodic boundary conditions when calculating angles
- [
False
] this is important when connecting vectors between atoms might require minimum image conventionUses cython implementation.
Changed in version 0.9.0: Added pbc option (default
False
)
-
dump_contents
()¶ Return a tuple of tuples which define the contents of this TopologyGroup in terms of the atom numbers, (0 based index within u.atoms)
This format should be identical to the original contents of the entries in universe._topology. Note that because bonds are sorted as they are initialised, the order that atoms are defined in each entry might be reversed.
New in version 0.9.0.
Changed in version 0.10.0: Renamed from “dump_contents” to “to_indices”
-
classmethod
from_indices
(bondlist, atomgroup, bondclass=None, guessed=True, remove_duplicates=False)[source]¶ Initiate from a list of indices.
Parameters: - *bondlist* – A list of lists of indices. For example [(0, 1), (1, 2)] Note that these indices refer to the index of the Atoms within the supplied AtomGroup, not their global index.
- *atomgroup* – An AtomGroup which the indices from bondlist will be used on.
Keywords: - bondclass
The Class of the topology object to be made. If missing this will try and be guessed according to the number of indices in each record.
- guessed
Whether or not the bonds were guessed. [
True
]- remove_duplicates
Sort through items and make sure that no duplicates are created [
False
]
New in version 0.10.0.
-
selectBonds
(selection)¶ Retrieves a selection from this topology group based on types.
-
to_indices
()[source]¶ Return a tuple of tuples which define the contents of this TopologyGroup in terms of the atom numbers, (0 based index within u.atoms)
This format should be identical to the original contents of the entries in universe._topology. Note that because bonds are sorted as they are initialised, the order that atoms are defined in each entry might be reversed.
New in version 0.9.0.
Changed in version 0.10.0: Renamed from “dump_contents” to “to_indices”
-
topDict
¶ Returns the TopologyDict for this topology group.
This is used for the select_bonds method when fetching a certain type of bond.
This is a cached property so will be generated the first time it is accessed.
-
-
class
MDAnalysis.core.topologyobjects.
TopologyObject
(atoms, is_guessed=False)[source]¶ Base class for all Topology items.
Defines the behaviour by which Bonds/Angles/etc in MDAnalysis should behave.
New in version 0.9.0.
Changed in version 0.10.0: All TopologyObject now keep track of if they were guessed or not via the
is_guessed
managed property.New in version 0.11.0: Added the value method to return the size of the object
-
__contains__
(other)[source]¶ Check whether an atom is in this
TopologyObject
-
indices
¶ Tuple of indices describing this object
New in version 0.10.0.
-
is_guessed
¶ True
if the bond was guessed.See also
guess_bonds()
guess_angles()
andguess_dihedrals()
-
type
¶ Type of the bond as a tuple
When comparing types, it is important to consider the reverse of the type too, i.e.:
a.type == b.type or a.type == b.type[::-1]
-