AtomGroup

A Universe contains all particles in the molecular system. MDAnalysis calls a particle an Atom, regardless of whether it really is (e.g. it may be a united-atom particle or coarse-grained bead). Atoms are grouped with an AtomGroup; the ‘master’ AtomGroup of a Universe is accessible at Universe.atoms.

Note

The AtomGroup is probably the most important object in MDAnalysis. Virtually everything can be accessed through an AtomGroup.

Creating an AtomGroup

Atom selection language

AtomGroup instances are typically created with Universe.select_atoms or by manipulating another AtomGroup, e.g. by slicing.

In [1]: import MDAnalysis as mda

In [2]: from MDAnalysis.tests.datafiles import PDB

In [3]: u = mda.Universe(PDB)

In [4]: u.select_atoms('resname ARG')
Out[4]: <AtomGroup with 312 atoms>

See Atom selection language for more information.

Indexing and slicing

An AtomGroup can be indexed and sliced like a list:

In [5]: print(u.atoms[0])
<Atom 1: N of type N of resname MET, resid 1 and segid SYSTEM and altLoc >

Slicing returns another AtomGroup. The below code returns an AtomGroup of every second element from the first to the 6th element, corresponding to indices 0, 2, and 4.

In [6]: ag = u.atoms[0:6:2]

In [7]: ag.indices
Out[7]: array([0, 2, 4])

MDAnalysis also supports fancy indexing: passing a ndarray or a list.

In [8]: indices = [0, 3, -1, 10, 3]

In [9]: u.atoms[indices].indices
Out[9]: array([    0,     3, 47680,    10,     3])

Boolean indexing allows you to pass in an array of True or False values to create a new AtomGroup from another. The array must be the same length as the original AtomGroup. This allows you to select atoms on conditions.

In [10]: arr = u.atoms.resnames == 'ARG'

In [11]: len(arr) == len(u.atoms)

In [12]: arr
Out[12]: Out[11]: array([False, False, False, ..., False, False, False])

In [13]: u.atoms[arr]

Group operators and set methods

MDAnalysis supports a number of ways to compare AtomGroups or construct a new one: group operators (e.g. concatenate(), subtract()) and set methods (e.g. union(), difference()). Group operators achieve a similar outcome to set methods. However, a key difference is that concatenate() and subtract() preserve the order of the atoms and any duplicates. union() and difference() return an AtomGroup where each atom is unique, and ordered by its topology index.

In [14]: ag1 = u.atoms[1:6]

In [15]: ag2 = u.atoms[8:3:-1]

In [16]: concat = ag1 + ag2

In [17]: concat.indices
Out[17]: array([1, 2, 3, 4, 5, 8, 7, 6, 5, 4])

In [18]: union = ag1 | ag2

In [19]: union.indices
Out[19]: array([1, 2, 3, 4, 5, 6, 7, 8])

Available operators

Unlike set methods and atom selection language, concatenation and subtraction keep the order of the atoms as well as duplicates.

Operation

Equivalent

Result

len(s)

number of atoms in the group

s == t

test if s and t contain the same elements in the same order

s.concatenate(t)

s + t

new Group with elements from s and from t

s.subtract(t)

new Group with elements from s that are not in t

Available set methods

Each of these methods create groups that are sorted sets of unique Atoms.

Operation

Equivalent

Result

s.isdisjoint(t)

True if s and t do not share elements

s.issubset(t)

test if all elements of s are part of t

s.is_strict_subset(t)

test if all elements of s are part of t, and s != t

s.issuperset(t)

test if all elements of t are part of s

s.is_strict_superset(t)

test if all elements of t are part of s, and s != t

s.union(t)

s | t

new Group with elements from both s and t

s.intersection(t)

s & t

new Group with elements common to s and t

s.difference(t)

s - t

new Group with elements of s that are not in t

s.symmetric_difference(t)

s ^ t

new Group with elements that are part of s or t but not both

Groupby and split

An AtomGroup can be constructed from another by separating atoms by properties.

AtomGroup.split can create a list of AtomGroups by splitting another AtomGroup by the ‘level’ of connectivity: one of atom, residue, molecule, or segment.

In [20]: ag1 = u.atoms[:100]

In [21]: ag1.split('residue')
Out[21]: 
[<AtomGroup with 19 atoms>,
 <AtomGroup with 24 atoms>,
 <AtomGroup with 19 atoms>,
 <AtomGroup with 19 atoms>,
 <AtomGroup with 19 atoms>]

An AtomGroup can also be separated according to values of topology attributes to produce a dictionary of {value:AtomGroup}.

In [22]: u.atoms.groupby('masses')
Out[22]: 
{32.06: <AtomGroup with 7 atoms>,
 1.008: <AtomGroup with 23853 atoms>,
 0.0: <AtomGroup with 11084 atoms>,
 12.011: <AtomGroup with 1040 atoms>,
 14.007: <AtomGroup with 289 atoms>,
 15.999: <AtomGroup with 11404 atoms>,
 22.98977: <AtomGroup with 4 atoms>}

Passing in multiple attributes groups them in order:

In [23]: u.select_atoms('resname SOL NA+').groupby(['masses', 'resnames'])
Out[23]: 
{(0.0, 'SOL'): <AtomGroup with 11084 atoms>,
 (1.008, 'SOL'): <AtomGroup with 22168 atoms>,
 (22.98977, 'NA+'): <AtomGroup with 4 atoms>,
 (15.999, 'SOL'): <AtomGroup with 11084 atoms>}

Constructing from Atoms

An AtomGroup can be created from an iterable of Atom instances:

In [24]: atom1 = u.atoms[4]

In [25]: atom2 = u.atoms[6]

In [26]: atom3 = u.atoms[2]

In [27]: ag = mda.AtomGroup([atom1, atom2, atom3])

In [28]: print(ag)
<AtomGroup [<Atom 5: CA of type C of resname MET, resid 1 and segid SYSTEM and altLoc >, <Atom 7: CB of type C of resname MET, resid 1 and segid SYSTEM and altLoc >, <Atom 3: H2 of type H of resname MET, resid 1 and segid SYSTEM and altLoc >]>

A neat shortcut for this is to simply add an Atom to another Atom or AtomGroup:

In [29]: ag = atom1 + atom2

In [30]: print(ag)
<AtomGroup [<Atom 5: CA of type C of resname MET, resid 1 and segid SYSTEM and altLoc >, <Atom 7: CB of type C of resname MET, resid 1 and segid SYSTEM and altLoc >]>

In [31]: ag += atom3

In [32]: print(ag)
<AtomGroup [<Atom 5: CA of type C of resname MET, resid 1 and segid SYSTEM and altLoc >, <Atom 7: CB of type C of resname MET, resid 1 and segid SYSTEM and altLoc >, <Atom 3: H2 of type H of resname MET, resid 1 and segid SYSTEM and altLoc >]>

An alternative method is to provide a list of indices and the Universe that the Atoms belong to:

In [33]: ag = mda.AtomGroup([4, 6, 2], u)

In [34]: print(ag)
<AtomGroup [<Atom 5: CA of type C of resname MET, resid 1 and segid SYSTEM and altLoc >, <Atom 7: CB of type C of resname MET, resid 1 and segid SYSTEM and altLoc >, <Atom 3: H2 of type H of resname MET, resid 1 and segid SYSTEM and altLoc >]>

Order and uniqueness

These methods of creating an AtomGroup result in a sorted, unique list of atoms:

These methods return a user-ordered AtomGroup that can contain duplicates:

Empty AtomGroups

MDAnalysis can also work with empty AtomGroups:

In [35]: null = u.atoms[[]]

In [36]: null
Out[36]: <AtomGroup with 0 atoms>

The above is the same as creating an AtomGroup from an empty list and a Universe.

In [37]: mda.AtomGroup([], u)
Out[37]: <AtomGroup with 0 atoms>

Each method of creating an AtomGroup can also be used to create an empty one. For example, using selection language:

In [38]: u.select_atoms("resname DOES_NOT_EXIST")
Out[38]: <AtomGroup with 0 atoms>

and indexing:

In [39]: u.atoms[6:6]
Out[39]: <AtomGroup with 0 atoms>

or set methods:

In [40]: u.atoms - u.atoms
Out[40]: <AtomGroup with 0 atoms>

Empty AtomGroups have a length of 0 and evaluate to False in a boolean context.

In [41]: bool(null)
Out[41]: False

This allows analysis methods to skip over empty AtomGroups instead of raising an error, which is helpful as occasionally empty AtomGroups can arise from selection logic that is too restrictive (e.g. geometric selections).

Dynamically updating AtomGroups

A normal AtomGroup is static, and the atoms within it do not change as the trajectory frame changes. Several methods require dynamically updating AtomGroups. These are typically created using atom selection language. See Dynamic selections for more information.

Methods

Most of the analysis functionality in MDAnalysis is implemented in the analysis module, but many interesting methods can be accessed from an AtomGroup directly. For example, Bonds, Angles, Dihedrals and ImproperDihedrals can be created from AtomGroups. Providing that required topology attributes are present, a number of analysis methods are also available to a AtomGroup, ResidueGroup, and SegmentGroup.