11.1.2. Core objects: Containers — MDAnalysis.core.groups

The Universe instance contains all the particles in the system (which MDAnalysis calls Atom). Groups of atoms are handled as AtomGroup instances. The AtomGroup is probably the most important object in MDAnalysis because virtually everything can be accessed through it. AtomGroup instances can be easily created (e.g., from an AtomGroup.select_atoms() selection or simply by slicing).

For convenience, chemically meaningful groups of Atoms such as a Residue or a Segment (typically a whole molecule or all of the solvent) also exist as containers, as well as groups of these units (ResidueGroup, SegmentGroup).

11.1.2.1. Classes

11.1.2.1.1. Collections

class MDAnalysis.core.groups.AtomGroup(*args)[source]

An ordered array of atoms.

Can be initiated from an iterable of Atoms:

ag = AtomGroup([Atom1, Atom2, Atom3])

Or from providing a list of indices and the Universe it should belong to:

ag = AtomGroup([72, 14, 25], u)

Alternatively, an AtomGroup is generated by indexing/slicing another AtomGroup, such as the group of all Atoms in the Universe at MDAnalysis.core.universe.Universe.atoms.

An AtomGroup can be indexed and sliced like a list:

ag[0], ag[-1]

will return the first and the last Atom in the group whereas the slice:

ag[0:6:2]

returns an AtomGroup of every second element, corresponding to indices 0, 2, and 4.

It also supports “advanced slicing” when the argument is a numpy.ndarray or a list:

aslice = [0, 3, -1, 10, 3]
ag[aslice]

will return a new AtomGroup of Atoms with those indices in the old AtomGroup.

Finally, AtomGroups can be created from a selection. See select_atoms().

Note

AtomGroups originating from a selection are sorted and duplicate elements are removed. This is not true for AtomGroups produced by slicing. Thus, slicing can be used when the order of atoms is crucial (for instance, in order to define angles or dihedrals).

AtomGroups can be compared and combined using group operators. For instance, AtomGroups can be concatenated using + or concatenate():

ag_concat = ag1 + ag2  # or ag_concat = ag1.concatenate(ag2)

When groups are concatenated, the order of the Atoms is conserved. If Atoms appear several times in one of the groups, the duplicates are kept in the resulting group. On the contrary to concatenate(), union() treats the AtomGroups as sets so that duplicates are removed from the resulting group, and Atoms are ordered. The | operator is synomymous to union():

ag_union = ag1 | ag2  # or ag_union = ag1.union(ag2)

The opposite operation to concatenate() is subtract(). This method creates a new group with all the Atoms of the group that are not in a given other group; the order of the Atoms is kept, and so are duplicates. difference() is the set version of subtract(). The resulting group is sorted and deduplicated.

All set methods are listed in the table below. These methods treat the groups as sorted and deduplicated sets of 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

The following methods keep the order of the atoms as well as duplicates.

Operation Equivalent Result
len(s)   number of elements (atoms, residues or segment) 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

The in operator allows to test if an Atom is in the AtomGroup.

AtomGroup instances are always bound to a MDAnalysis.core.universe.Universe. They cannot exist in isolation.

See also

MDAnalysis.core.universe.Universe

Instant selectors of AtomGroup will be removed in the 1.0 release.
Removed instant selectors, use select_atoms(‘name …’) to select atoms by name.
accumulate(attribute, function=<function sum>, compound='group')

Accumulates the attribute associated with (compounds of) the group.

Accumulates the attribute of Atoms in the group. The accumulation per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly. By default, the method sums up all attributes per compound, but any function that takes an array and returns an acuumulation over a given axis can be used. For multi-dimensional input arrays, the accumulation is performed along the first axis.

Parameters:
  • attribute (str or array_like) – Attribute or array of values to accumulate. If a numpy.ndarray (or compatible) is provided, its first dimension must have the same length as the total number of atoms in the group.
  • function (callable, optional) – The function performing the accumulation. It must take the array of attribute values to accumulate as its only positional argument and accept an (optional) keyword argument axis allowing to specify the axis along which the accumulation is performed.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the accumulation of all attributes associated with atoms in the group will be returned as a single value. Otherwise, the accumulation of the attributes per Segment, Residue, molecule, or fragment will be returned as a 1d array. Note that, in any case, only the Atoms belonging to the group will be taken into account.
Returns:

Acuumulation of the attribute. If compound is set to 'group', the first dimension of the attribute array will be contracted to a single value. If compound is set to 'segments', 'residues', 'molecules', or 'fragments', the length of the first dimension will correspond to the number of compounds. In all cases, the other dimensions of the returned array will be of the original shape (without the first dimension).

Return type:

float or numpy.ndarray

Raises:
  • ValueError – If the length of a provided attribute array does not correspond to the number of atoms in the group.
  • ValueError – If compound is not one of 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums), or if compound is 'fragments' but the topology doesn’t contain bonds.

Examples

To find the total charge of a given AtomGroup:

>>> sel = u.select_atoms('prop mass > 4.0')
>>> sel.accumulate('charges')

To find the total mass per residue of all CA Atoms:

>>> sel = u.select_atoms('name CA')
>>> sel.accumulate('masses', compound='residues')

To find the maximum atomic charge per fragment of a given AtomGroup:

>>> sel.accumulate('charges', compound="fragments", function=np.max)

New in version 0.20.0.

angle

This AtomGroup represented as an MDAnalysis.core.topologyobjects.Angle object

Raises:ValueError – If the AtomGroup is not length 3

New in version 0.11.0.

atoms

The AtomGroup itself.

See also

copy
return a true copy of the AtomGroup

Changed in version 0.19.0: In previous versions, this returned a copy, but now the AtomGroup itself is returned. This should not affect any code but only speed up calculations.

bbox(pbc=False)

Return the bounding box of the selection.

The lengths A,B,C of the orthorhombic enclosing box are

L = AtomGroup.bbox()
A,B,C = L[1] - L[0]
Parameters:pbc (bool, optional) – If True, move all Atoms to the primary unit cell before calculation. [False]
Returns:corners – 2x3 array giving corners of bounding box as [[xmin, ymin, zmin], [xmax, ymax, zmax]].
Return type:numpy.ndarray

New in version 0.7.2.

Changed in version 0.8: Added pbc keyword

Changed in version 1.0.0: Removed flags affecting default behaviour

bond

This AtomGroup represented as a MDAnalysis.core.topologyobjects.Bond object

Raises:ValueError – If the AtomGroup is not length 2

New in version 0.11.0.

bsphere(pbc=False)

Return the bounding sphere of the selection.

The sphere is calculated relative to the center of geometry.

Parameters:pbc (bool, optional) – If True, move all atoms to the primary unit cell before calculation. [False]
Returns:
  • R (float) – Radius of the bounding sphere.
  • center (numpy.ndarray) – Coordinates of the sphere center as [xcen, ycen, zcen].

New in version 0.7.3.

Changed in version 0.8: Added pbc keyword

center(weights, pbc=False, compound='group', unwrap=False)

Weighted center of (compounds of) the group

Computes the weighted center of Atoms in the group. Weighted centers per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly. If the weights of a compound sum up to zero, the coordinates of that compound’s weighted center will be nan (not a number).

Parameters:
  • weights (array_like or None) – Weights to be used. Setting weights=None is equivalent to passing identical weights for all atoms of the group.
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments', 'residues', 'molecules', or 'fragments', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default [False].
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the weighted center of all atoms in the group will be returned as a single position vector. Else, the weighted centers of each Segment, Residue, molecule, or fragment will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the weighted center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments', 'residues', 'molecules', or 'fragments', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Raises:
  • ValueError – If compound is not one of 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • ValueError – If both ‘pbc’ and ‘unwrap’ set to true.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums) or if compound is 'fragments' but the topology doesn’t contain bonds.

Examples

To find the center of charge of a given AtomGroup:

>>> sel = u.select_atoms('prop mass > 4.0')
>>> sel.center(sel.charges)

To find the centers of mass per residue of all CA Atoms:

>>> sel = u.select_atoms('name CA')
>>> sel.center(sel.masses, compound='residues')

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

center_of_geometry(pbc=False, compound='group', unwrap=False)

Center of geometry of (compounds of) the group.

Computes the center of geometry (a.k.a. centroid) of Atoms in the group. Centers of geometry per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly.

Parameters:
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments' or 'residues', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default False.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the center of geometry of all Atoms in the group will be returned as a single position vector. Else, the centers of geometry of each Segment or Residue will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the geometric center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments' or 'residues', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Changed in version 0.8: Added pbc keyword

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

centroid(pbc=False, compound='group', unwrap=False)

Center of geometry of (compounds of) the group.

Computes the center of geometry (a.k.a. centroid) of Atoms in the group. Centers of geometry per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly.

Parameters:
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments' or 'residues', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default False.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the center of geometry of all Atoms in the group will be returned as a single position vector. Else, the centers of geometry of each Segment or Residue will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the geometric center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments' or 'residues', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Changed in version 0.8: Added pbc keyword

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

cmap

This AtomGroup represented as an MDAnalysis.core.topologyobjects.CMap object

Raises:ValueError – If the AtomGroup is not length 5

New in version 1.0.0.

concatenate(other)

Concatenate with another Group or Component of the same level.

Duplicate entries and original order is preserved. It is synomymous to the + operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with elements of self and other concatenated
Return type:Group

Example

The order of the original contents (including duplicates) are preserved when performing a concatenation.

>>> ag1 = u.select_atoms('name O')
>>> ag2 = u.select_atoms('name N')
>>> ag3 = ag1 + ag2  # or ag1.concatenate(ag2)
>>> ag3[:3].names
array(['O', 'O', 'O'], dtype=object)
>>> ag3[-3:].names
array(['N', 'N', 'N'], dtype=object)

New in version 0.16.0.

convert_to(package)[source]

Convert AtomGroup to a structure from another Python package.

Example

The code below converts a Universe to a parmed.structure.Structure.

>>> import MDAnalysis as mda
>>> from MDAnalysis.tests.datafiles import GRO
>>> u = mda.Universe(GRO)
>>> parmed_structure = u.atoms.convert_to('PARMED')
>>> parmed_structure
<Structure 47681 atoms; 11302 residues; 0 bonds; PBC (triclinic); NOT parametrized>
Parameters:package (str) – The name of the package to convert to, e.g. "PARMED"
Returns:An instance of the structure type from another package.
Return type:output
Raises:TypeError: – No converter was found for the required package

New in version 1.0.0.

copy()

Get another group identical to this one.

New in version 0.19.0.

difference(other)

Elements from this Group that do not appear in another

This method removes duplicate elements and sorts the result. As such, it is different from subtract(). difference() is synomymous to the - operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements of self that are not in other, without duplicate elements
Return type:Group

New in version 0.16.

dihedral

This AtomGroup represented as a Dihedral object

Raises:ValueError – If the AtomGroup is not length 4

New in version 0.11.0.

dimensions

Obtain a copy of the dimensions of the currently loaded Timestep

forces

Forces on each Atom in the AtomGroup.

A numpy.ndarray with shape=(n_atoms, 3) and dtype=numpy.float32.

The forces can be changed by assigning an array of the appropriate shape, i.e. either (n_atoms, 3) to assign individual forces or (3,) to assign the same force to all Atoms (e.g. ag.forces = array([0,0,0]) will give all Atoms a zero force).

Raises:NoDataError – If the Timestep does not contain forces.
groupby(topattrs)

Group together items in this group according to values of topattr

Parameters:topattrs (str or list) – One or more topology attributes to group components by. Single arguments are passed as a string. Multiple arguments are passed as a list.
Returns:Unique values of the multiple combinations of topology attributes as keys, Groups as values.
Return type:dict

Example

To group atoms with the same mass together:

>>> ag.groupby('masses')
{12.010999999999999: <AtomGroup with 462 atoms>,
 14.007: <AtomGroup with 116 atoms>,
 15.999000000000001: <AtomGroup with 134 atoms>}

To group atoms with the same residue name and mass together:

>>> ag.groupby(['resnames', 'masses'])
{('ALA', 1.008): <AtomGroup with 95 atoms>,
 ('ALA', 12.011): <AtomGroup with 57 atoms>,
 ('ALA', 14.007): <AtomGroup with 19 atoms>,
 ('ALA', 15.999): <AtomGroup with 19 atoms>},
 ('ARG', 1.008): <AtomGroup with 169 atoms>,
 ...}
>>> ag.groupby(['resnames', 'masses'])('ALA', 15.999)
 <AtomGroup with 19 atoms>

New in version 0.16.0.

Changed in version 0.18.0: The function accepts multiple attributes

guess_bonds(vdwradii=None)[source]

Guess bonds that exist within this AtomGroup and add them to the underlying universe.

Parameters:vdwradii (dict, optional) – Dict relating atom types: vdw radii

New in version 0.10.0.

Changed in version 0.20.2: Now applies periodic boundary conditions when guessing bonds.

improper

This AtomGroup represented as an MDAnalysis.core.topologyobjects.ImproperDihedral object

Raises:ValueError – If the AtomGroup is not length 4

New in version 0.11.0.

intersection(other)

Group of elements which are in both this Group and another

This method removes duplicate elements and sorts the result. It is synomymous to the & operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the common elements of self and other, without duplicate elements
Return type:Group

Example

Intersections can be used when the select atoms string would become too complicated. For example to find the water atoms which are within 4.0A of two segments:

>>> shell1 = u.select_atoms('resname SOL and around 4.0 segid 1')
>>> shell2 = u.select_atoms('resname SOL and around 4.0 segid 2')
>>> common = shell1 & shell2  # or shell1.intersection(shell2)

See also

union()

New in version 0.16.

is_strict_subset(other)

If this Group is a subset of another Group but not identical

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a strict subset of the other one
Return type:bool

New in version 0.16.

is_strict_superset(other)

If this Group is a superset of another Group but not identical

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a strict superset of the other one
Return type:bool

New in version 0.16.

isdisjoint(other)

If the Group has no elements in common with the other Group

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if the two Groups do not have common elements
Return type:bool

New in version 0.16.

issubset(other)

If all elements of this Group are part of another Group

Note that an empty group is a subset of any group of the same level.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a subset of the other one
Return type:bool

New in version 0.16.

issuperset(other)

If all elements of another Group are part of this Group

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a subset of the other one
Return type:bool

New in version 0.16.

isunique

Boolean indicating whether all components of the group are unique, i.e., the group contains no duplicates.

Examples

>>> ag = u.atoms[[2, 1, 2, 2, 1, 0]]
>>> ag
<AtomGroup with 6 atoms>
>>> ag.isunique
False
>>> ag2 = ag.unique
>>> ag2
<AtomGroup with 3 atoms>
>>> ag2.isunique
True

New in version 0.19.0.

ix

Unique indices of the components in the Group.

ix_array

Unique indices of the components in the Group.

For a Group, ix_array is the same as ix. This method gives a consistent API between components and groups.

See also

ix

n_atoms

Number of atoms in the AtomGroup.

Equivalent to len(self).

n_residues

Number of unique Residues present in the AtomGroup.

Equivalent to len(self.residues).

n_segments

Number of unique segments present in the AtomGroup.

Equivalent to len(self.segments).

pack_into_box(box=None, inplace=True)

Shift all Atoms in this group to the primary unit cell.

Parameters:
  • box (array_like) – Box dimensions, can be either orthogonal or triclinic information. Cell dimensions must be in an identical to format to those returned by MDAnalysis.coordinates.base.Timestep.dimensions, [lx, ly, lz, alpha, beta, gamma]. If None, uses these timestep dimensions.
  • inplace (bool) – True to change coordinates in place.
Returns:

coords – Shifted atom coordinates.

Return type:

numpy.ndarray

Notes

All atoms will be moved so that they lie between 0 and boxlength \(L_i\) in all dimensions, i.e. the lower left corner of the simulation box is taken to be at (0,0,0):

\[x_i' = x_i - \left\lfloor\frac{x_i}{L_i}\right\rfloor\]

The default is to take unit cell information from the underlying Timestep instance. The optional argument box can be used to provide alternative unit cell information (in the MDAnalysis standard format [Lx, Ly, Lz, alpha, beta, gamma]).

Works with either orthogonal or triclinic box types.

Note

pack_into_box() is identical to wrap() with all default keywords.

New in version 0.8.

positions

Coordinates of the Atoms in the AtomGroup.

A numpy.ndarray with shape=(n_atoms, 3) and dtype=numpy.float32.

The positions can be changed by assigning an array of the appropriate shape, i.e., either (n_atoms, 3) to assign individual coordinates, or (3,) to assign the same coordinate to all Atoms (e.g., ag.positions = array([0,0,0]) will move all Atoms to the origin).

Note

Changing positions is not reflected in any files; reading any frame from the trajectory will replace the change with that from the file except if the trajectory is held in memory, e.g., when the transfer_to_memory() method was used.

Raises:NoDataError – If the underlying Timestep does not contain positions.
residues

A sorted ResidueGroup of the unique Residues present in the AtomGroup.

rotate(R, point=(0, 0, 0))

Apply a rotation matrix R to the selection’s coordinates. \(\mathsf{R}\) is a 3x3 orthogonal matrix that transforms a vector \(\mathbf{x} \rightarrow \mathbf{x}'\):

\[\mathbf{x}' = \mathsf{R}\mathbf{x}\]

Atom coordinates are rotated in-place.

Parameters:
  • R (array_like) – 3x3 rotation matrix
  • point (array_like, optional) – Center of rotation
Returns:

Return type:

self

Notes

By default, rotates about the origin point=(0, 0, 0). To rotate a group g around its center of geometry, use g.rotate(R, point=g.center_of_geometry()).

See also

rotateby()
rotate around given axis and angle
MDAnalysis.lib.transformations()
module of all coordinate transforms
rotateby(angle, axis, point=None)

Apply a rotation to the selection’s coordinates.

Parameters:
  • angle (float) – Rotation angle in degrees.
  • axis (array_like) – Rotation axis vector.
  • point (array_like, optional) – Center of rotation. If None then the center of geometry of this group is used.
Returns:

Return type:

self

Notes

The transformation from current coordinates \(\mathbf{x}\) to new coordinates \(\mathbf{x}'\) is

\[\mathbf{x}' = \mathsf{R}\,(\mathbf{x}-\mathbf{p}) + \mathbf{p}\]

where \(\mathsf{R}\) is the rotation by angle around the axis going through point \(\mathbf{p}\).

See also

MDAnalysis.lib.transformations.rotation_matrix()

calculate()
math:mathsf{R}
segments

A sorted SegmentGroup of the unique segments present in the AtomGroup.

select_atoms(sel, *othersel, **selgroups)[source]

Select atoms from within this Group using a selection string.

Returns an AtomGroup sorted according to their index in the topology (this is to ensure that there are no duplicates, which can happen with complicated selections).

Parameters:
  • sel (str) – string of the selection, eg “name Ca”, see below for possibilities.
  • othersel (iterable of str) – further selections to perform. The results of these selections will be appended onto the results of the first.
  • periodic (bool (optional)) – for geometric selections, whether to account for atoms in different periodic images when searching
  • updating (bool (optional)) – force the selection to be re evaluated each time the Timestep of the trajectory is changed. See section on Dynamic selections below. [True]
  • **selgroups (keyword arguments of str: AtomGroup (optional)) – when using the “group” keyword in selections, groups are defined by passing them as keyword arguments. See section on preexisting selections below.
Raises:

TypeError – If the arbitrary groups passed are not of type MDAnalysis.core.groups.AtomGroup

Examples

All simple selection listed below support multiple arguments which are implicitly combined with an or operator. For example

>>> sel = universe.select_atoms('resname MET GLY')

is equivalent to

>>> sel = universe.select_atoms('resname MET or resname GLY')

Will select all atoms with a residue name of either MET or GLY.

Subselections can be grouped with parentheses.

>>> sel = universe.select_atoms("segid DMPC and not ( name H* O* )")
>>> sel
<AtomGroup with 3420 atoms>

Existing AtomGroup objects can be passed as named arguments, which will then be available to the selection parser.

>>> universe.select_atoms("around 10 group notHO", notHO=sel)
<AtomGroup with 1250 atoms>

Selections can be set to update automatically on frame change, by setting the updating keyword argument to True. This will return a UpdatingAtomGroup which can represent the solvation shell around another object.

>>> universe.select_atoms("resname SOL and around 2.0 protein", updating=True)
<Updating AtomGroup with 100 atoms>

Notes

If exact ordering of atoms is required (for instance, for angle() or dihedral() calculations) then one supplies selections separately in the required order. Also, when multiple AtomGroup instances are concatenated with the + operator, then the order of Atom instances is preserved and duplicates are not removed.

Selection syntax

The selection parser understands the following CASE SENSITIVE keywords:

Simple selections

protein, backbone, nucleic, nucleicbackbone
selects all atoms that belong to a standard set of residues; a protein is identfied by a hard-coded set of residue names so it may not work for esoteric residues.
segid seg-name
select by segid (as given in the topology), e.g. segid 4AKE or segid DMPC
resid residue-number-range
resid can take a single residue number or a range of numbers. A range consists of two numbers separated by a colon (inclusive) such as resid 1:5. A residue number (“resid”) is taken directly from the topology. If icodes are present in the topology, then these will be taken into account. Ie ‘resid 163B’ will only select resid 163 with icode B while ‘resid 163’ will select only residue 163. Range selections will also respect icodes, so ‘resid 162-163B’ will select all residues in 162 and those in 163 up to icode B.
resnum resnum-number-range
resnum is the canonical residue number; typically it is set to the residue id in the original PDB structure.
resname residue-name
select by residue name, e.g. resname LYS
name atom-name
select by atom name (as given in the topology). Often, this is force field dependent. Example: name CA (for C&alpha; atoms) or name OW (for SPC water oxygen)
type atom-type
select by atom type; this is either a string or a number and depends on the force field; it is read from the topology file (e.g. the CHARMM PSF file contains numeric atom types). It has non-sensical values when a PDB or GRO file is used as a topology
atom seg-name residue-number atom-name
a selector for a single atom consisting of segid resid atomname, e.g. DMPC 1 C2 selects the C2 carbon of the first residue of the DMPC segment
altloc alternative-location
a selection for atoms where alternative locations are available, which is often the case with high-resolution crystal structures e.g. resid 4 and resname ALA and altloc B selects only the atoms of ALA-4 that have an altloc B record.
moltype molecule-type
select by molecule type, e.g. moltype Protein_A. At the moment, only the TPR format defines the molecule type.
record_type record_type
for selecting either ATOM or HETATM from PDB-like files. e.g. select_atoms('name CA and not record_type HETATM')

Boolean

not
all atoms not in the selection, e.g. not protein selects all atoms that aren’t part of a protein
and, or
combine two selections according to the rules of boolean algebra, e.g. protein and not resname ALA LYS selects all atoms that belong to a protein, but are not in a lysine or alanine residue

Geometric

around distance selection
selects all atoms a certain cutoff away from another selection, e.g. around 3.5 protein selects all atoms not belonging to protein that are within 3.5 Angstroms from the protein
point x y z distance
selects all atoms within a cutoff of a point in space, make sure coordinate is separated by spaces, e.g. point 5.0 5.0 5.0  3.5 selects all atoms within 3.5 Angstroms of the coordinate (5.0, 5.0, 5.0)
prop [abs] property operator value
selects atoms based on position, using property x, y, or z coordinate. Supports the abs keyword (for absolute value) and the following operators: <, >, <=, >=, ==, !=. For example, prop z >= 5.0 selects all atoms with z coordinate greater than 5.0; prop abs z <= 5.0 selects all atoms within -5.0 <= z <= 5.0.
sphzone radius selection
Selects all atoms that are within radius of the center of geometry of selection
sphlayer inner radius outer radius selection
Similar to sphzone, but also excludes atoms that are within inner radius of the selection COG
cyzone externalRadius zMax zMin selection
selects all atoms within a cylindric zone centered in the center of geometry (COG) of a given selection, e.g. cyzone 15 4 -8 protein and resid 42 selects the center of geometry of protein and resid 42, and creates a cylinder of external radius 15 centered on the COG. In z, the cylinder extends from 4 above the COG to 8 below. Positive values for zMin, or negative ones for zMax, are allowed.
cylayer innerRadius externalRadius zMax zMin selection
selects all atoms within a cylindric layer centered in the center of geometry (COG) of a given selection, e.g. cylayer 5 10 10 -8 protein selects the center of geometry of protein, and creates a cylindrical layer of inner radius 5, external radius 10 centered on the COG. In z, the cylinder extends from 10 above the COG to 8 below. Positive values for zMin, or negative ones for zMax, are allowed.

Connectivity

byres selection
selects all atoms that are in the same segment and residue as selection, e.g. specify the subselection after the byres keyword
bonded selection
selects all atoms that are bonded to selection eg: select name H and bonded name O selects only hydrogens bonded to oxygens

Index

bynum index-range
selects all atoms within a range of (1-based) inclusive indices, e.g. bynum 1 selects the first atom in the universe; bynum 5:10 selects atoms 5 through 10 inclusive. All atoms in the Universe are consecutively numbered, and the index runs from 1 up to the total number of atoms.
index index-range
selects all atoms within a range of (0-based) inclusive indices, e.g. index 0 selects the first atom in the universe; index 5:10 selects atoms 6 through 11 inclusive. All atoms in the Universe are consecutively numbered, and the index runs from 0 up to the total number of atoms - 1.

Preexisting selections

group group-name
selects the atoms in the AtomGroup passed to the function as a keyword argument named group-name. Only the atoms common to group-name and the instance select_atoms() was called from will be considered, unless group is preceded by the global keyword. group-name will be included in the parsing just by comparison of atom indices. This means that it is up to the user to make sure the group-name group was defined in an appropriate Universe.
global selection
by default, when issuing select_atoms() from an AtomGroup, selections and subselections are returned intersected with the atoms of that instance. Prefixing a selection term with global causes its selection to be returned in its entirety. As an example, the global keyword allows for lipids.select_atoms("around 10 global protein") — where lipids is a group that does not contain any proteins. Were global absent, the result would be an empty selection since the protein subselection would itself be empty. When issuing select_atoms() from a Universe, global is ignored.
Dynamic selections
If select_atoms() is invoked with named argument updating set to True, an UpdatingAtomGroup instance will be returned, instead of a regular AtomGroup. It behaves just like the latter, with the difference that the selection expressions are re-evaluated every time the trajectory frame changes (this happens lazily, only when the UpdatingAtomGroup is accessed so that there is no redundant updating going on). Issuing an updating selection from an already updating group will cause later updates to also reflect the updating of the base group. A non-updating selection or a slicing operation made on an UpdatingAtomGroup will return a static AtomGroup, which will no longer update across frames.

Changed in version 0.7.4: Added resnum selection.

Changed in version 0.8.1: Added group and fullgroup selections.

Changed in version 0.13.0: Added bonded selection.

Changed in version 0.16.0: Resid selection now takes icodes into account where present.

Changed in version 0.16.0: Updating selections now possible by setting the updating argument.

Changed in version 0.17.0: Added moltype and molnum selections.

Changed in version 0.19.0: Added strict type checking for passed groups. Added periodic kwarg (default True)

Changed in version 0.19.2: Empty sel string now returns an empty Atom group.

Changed in version 1.0.0: The fullgroup selection has now been removed in favor of the equivalent global group selection. Removed flags affecting default behaviour for periodic selections; periodic are now on by default (as with default flags)

split(level)[source]

Split AtomGroup into a list of AtomGroups by level.

Parameters:level ({'atom', 'residue', 'molecule', 'segment'}) –

New in version 0.9.0.

Changed in version 0.17.0: Added the ‘molecule’ level.

subtract(other)

Group with elements from this Group that don’t appear in other

The original order of this group is kept, as well as any duplicate elements. If an element of this Group is duplicated and appears in the other Group or Component, then all the occurences of that element are removed from the returned Group.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements of self that are not in other, conserves order and duplicates.
Return type:Group

Example

Unlike difference() this method will not sort or remove duplicates.

>>> ag1 = u.atoms[[3, 3, 2, 2, 1, 1]]
>>> ag2 = u.atoms[2]
>>> ag3 = ag1 - ag2  # or ag1.subtract(ag2)
>>> ag1.indices
array([3, 3, 1, 1])

New in version 0.16.

symmetric_difference(other)

Group of elements which are only in one of this Group or another

This method removes duplicate elements and the result is sorted. It is synomym to the ^ operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements that are in self or in other but not in both, without duplicate elements
Return type:Group

Example

>>> ag1 = u.atoms[[0, 1, 5, 3, 3, 2]]
>>> ag2 = u.atoms[[4, 4, 6, 2, 3, 5]]
>>> ag3 = ag1 ^ ag2  # or ag1.symmetric_difference(ag2)
>>> ag3.indices  # 0 and 1 are only in ag1, 4 and 6 are only in ag2
[0, 1, 4, 6]

See also

difference()

New in version 0.16.

transform(M)

Apply homogenous transformation matrix M to the coordinates.

Atom coordinates are rotated and translated in-place.

Parameters:M (array_like) – 4x4 matrix with the rotation in R = M[:3, :3] and the translation in t = M[:3, 3].
Returns:
Return type:self

See also

MDAnalysis.lib.transformations()
module of all coordinate transforms

Notes

The rotation \(\mathsf{R}\) is about the origin and is applied before the translation \(\mathbf{t}\):

\[\mathbf{x}' = \mathsf{R}\mathbf{x} + \mathbf{t}\]
translate(t)

Apply translation vector t to the selection’s coordinates.

Atom coordinates are translated in-place.

Parameters:t (array_like) – vector to translate coordinates with
Returns:
Return type:self

See also

MDAnalysis.lib.transformations()
module of all coordinate transforms

Notes

The method applies a translation to the AtomGroup from current coordinates \(\mathbf{x}\) to new coordinates \(\mathbf{x}'\):

\[\mathbf{x}' = \mathbf{x} + \mathbf{t}\]
ts

Temporary Timestep that contains the selection coordinates.

A Timestep instance, which can be passed to a trajectory writer.

If ts is modified then these modifications will be present until the frame number changes (which typically happens when the underlying trajectory frame changes).

It is not possible to assign a new Timestep to the AtomGroup.ts attribute; change attributes of the object.

union(other)

Group of elements either in this Group or another

On the contrary to concatenation, this method sort the elements and removes duplicate ones. It is synomymous to the | operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the combined elements of self and other, without duplicate elements
Return type:Group

Example

In contrast to concatenate(), any duplicates are dropped and the result is sorted.

>>> ag1 = u.select_atoms('name O')
>>> ag2 = u.select_atoms('name N')
>>> ag3 = ag1 | ag2  # or ag1.union(ag2)
>>> ag3[:3].names
array(['N', 'O', 'N'], dtype=object)

New in version 0.16.

unique

An AtomGroup containing sorted and unique Atoms only.

If the AtomGroup is unique, this is the group itself.

Examples

>>> ag = u.atoms[[2, 1, 2, 2, 1, 0]]
>>> ag
<AtomGroup with 6 atoms>
>>> ag.ix
array([2, 1, 2, 2, 1, 0])
>>> ag2 = ag.unique
>>> ag2
<AtomGroup with 3 atoms>
>>> ag2.ix
array([0, 1, 2])
>>> ag2.unique is ag2
True

New in version 0.16.0.

Changed in version 0.19.0: If the AtomGroup is already unique, AtomGroup.unique now returns the group itself instead of a copy.

universe

The underlying Universe the group belongs to.

unwrap(compound='fragments', reference='com', inplace=True)

Move atoms of this group so that bonds within the group’s compounds aren’t split across periodic boundaries.

This function is most useful when atoms have been packed into the primary unit cell, causing breaks mid-molecule, with the molecule then appearing on either side of the unit cell. This is problematic for operations such as calculating the center of mass of the molecule.

+-----------+       +-----------+
|           |       |           |
| 6       3 |       |         3 | 6
| !       ! |       |         ! | !
|-5-8   1-2-|  ==>  |       1-2-|-5-8
| !       ! |       |         ! | !
| 7       4 |       |         4 | 7
|           |       |           |
+-----------+       +-----------+
Parameters:
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – Which type of component to unwrap. Note that, in any case, all atoms within each compound must be interconnected by bonds, i.e., compounds must correspond to (parts of) molecules.
  • reference ({'com', 'cog', None}, optional) – If 'com' (center of mass) or 'cog' (center of geometry), the unwrapped compounds will be shifted so that their individual reference point lies within the primary unit cell. If None, no such shift is performed.
  • inplace (bool, optional) – If True, coordinates are modified in place.
Returns:

coords – Unwrapped atom coordinate array of shape (n, 3).

Return type:

numpy.ndarray

Raises:
  • NoDataError – If compound is 'molecules' but the underlying topology does not contain molecule information, or if reference is 'com' but the topology does not contain masses.
  • ValueError – If reference is not one of 'com', 'cog', or None, or if reference is 'com' and the total mass of any compound is zero.

Note

Be aware of the fact that only atoms belonging to the group will be unwrapped! If you want entire molecules to be unwrapped, make sure that all atoms of these molecules are part of the group.

An AtomGroup containing all atoms of all fragments in the group ag can be created with:

all_frag_atoms = sum(ag.fragments)

See also

make_whole(), wrap(), pack_into_box(), apply_PBC()

New in version 0.20.0.

ureybradley

This AtomGroup represented as an MDAnalysis.core.topologyobjects.UreyBradley object

Raises:ValueError – If the AtomGroup is not length 2

New in version 1.0.0.

velocities

Velocities of the Atoms in the AtomGroup.

A numpy.ndarray with shape=(n_atoms, 3) and dtype=numpy.float32.

The velocities can be changed by assigning an array of the appropriate shape, i.e. either (n_atoms, 3) to assign individual velocities or (3,) to assign the same velocity to all Atoms (e.g. ag.velocities = array([0,0,0]) will give all Atoms zero velocity).

Raises:NoDataError – If the underlying Timestep does not contain velocities.
wrap(compound='atoms', center='com', box=None, inplace=True)

Shift the contents of this group back into the primary unit cell according to periodic boundary conditions.

Specifying a compound will keep the Atoms in each compound together during the process. If compound is different from 'atoms', each compound as a whole will be shifted so that its center lies within the primary unit cell.

Parameters:
  • compound ({'atoms', 'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – Which type of component to keep together during wrapping. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • center ({'com', 'cog'}) – How to define the center of a given group of atoms. If compound is 'atoms', this parameter is meaningless and therefore ignored.
  • box (array_like, optional) –

    The unitcell dimensions of the system, which can be orthogonal or triclinic and must be provided in the same format as returned by MDAnalysis.coordinates.base.Timestep.dimensions:

    [lx, ly, lz, alpha, beta, gamma].

    If None, uses the dimensions of the current time step.

  • inplace (bool, optional) – If True, coordinates will be changed in place.
Returns:

Array of wrapped atom coordinates of dtype np.float32 and shape (len(self.atoms.n_atoms), 3)

Return type:

numpy.ndarray

Raises:
  • ValueError – If compound is not one of 'atoms', 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums) or if compound is 'fragments' but the topology doesn’t contain bonds or if center is 'com' but the topology doesn’t contain masses.

Notes

All atoms of the group will be moved so that the centers of its compounds lie within the primary periodic image. For orthorhombic unit cells, the primary periodic image is defined as the half-open interval \([0,L_i)\) between \(0\) and boxlength \(L_i\) in all dimensions \(i\in\{x,y,z\}\), i.e., the origin of the of the simulation box is taken to be at the origin \((0,0,0)\) of the euclidian coordinate system. A compound center residing at position \(x_i\) in dimension \(i\) will be shifted to \(x_i'\) according to

\[x_i' = x_i - \left\lfloor\frac{x_i}{L_i}\right\rfloor\,.\]

When specifying a compound, the translation is calculated based on each compound. The same translation is applied to all atoms within this compound, meaning it will not be broken by the shift. This might however mean that not all atoms of a compound will be inside the unit cell after wrapping, but rather will be the center of the compound.

Be aware of the fact that only atoms belonging to the group will be taken into account!

center allows to define how the center of each group is computed. This can be either 'com' for center of mass, or 'cog' for center of geometry.

box allows a unit cell to be given for the transformation. If not specified, the dimensions information from the current Timestep will be used.

Note

AtomGroup.wrap() is currently faster than ResidueGroup.wrap() or SegmentGroup.wrap().

See also

pack_into_box(), unwrap(), MDanalysis.lib.distances.apply_PBC()

New in version 0.9.2.

Changed in version 0.20.0: The method only acts on atoms belonging to the group and returns the wrapped positions as a numpy.ndarray. Added optional argument inplace.

write(filename=None, file_format=None, filenamefmt='{trjname}_{frame}', frames=None, **kwargs)[source]

Write AtomGroup to a file.

The output can either be a coordinate file or a selection, depending on the format.

Examples

>>> ag = u.atoms
>>> ag.write('selection.ndx')  # Write a gromacs index file
>>> ag.write('coordinates.pdb')  # Write the current frame as PDB
>>> # Write the trajectory in XTC format
>>> ag.write('trajectory.xtc', frames='all')
>>> # Write every other frame of the trajectory in PBD format
>>> ag.write('trajectory.pdb', frames=u.trajectory[::2])
Parameters:
  • filename (str, optional) – None: create TRJNAME_FRAME.FORMAT from filenamefmt [None]
  • file_format (str, optional) – The name or extension of a coordinate, trajectory, or selection file format such as PDB, CRD, GRO, VMD (tcl), PyMol (pml), Gromacs (ndx) CHARMM (str) or Jmol (spt); case-insensitive [PDB]
  • filenamefmt (str, optional) – format string for default filename; use substitution tokens ‘trjname’ and ‘frame’ [“%(trjname)s_%(frame)d”]
  • bonds (str, optional) – how to handle bond information, especially relevant for PDBs. "conect": write only the CONECT records defined in the original file. "all": write out all bonds, both the original defined and those guessed by MDAnalysis. None: do not write out bonds. Default is "conect".
  • frames (array-like or slice or FrameIteratorBase or str, optional) – An ensemble of frames to write. The ensemble can be an list or array of frame indices, a mask of booleans, an instance of slice, or the value returned when a trajectory is indexed. By default, frames is set to None and only the current frame is written. If frames is set to “all”, then all the frame from trajectory are written.

Changed in version 0.9.0: Merged with write_selection. This method can now write both selections out.

Changed in version 0.19.0: Can write multiframe trajectories with the ‘frames’ argument.

class MDAnalysis.core.groups.ResidueGroup(*args)[source]

ResidueGroup base class.

This class is used by a Universe for generating its Topology-specific ResidueGroup class. All the TopologyAttr components are obtained from GroupBase, so this class only includes ad-hoc methods specific to ResidueGroups.

ResidueGroups can be compared and combined using group operators. See the list of these operators on GroupBase.

Deprecated since version 0.16.2: Instant selectors of Segments will be removed in the 1.0 release.

Changed in version 1.0.0: Removed instant selectors, use select_atoms instead

accumulate(attribute, function=<function sum>, compound='group')

Accumulates the attribute associated with (compounds of) the group.

Accumulates the attribute of Atoms in the group. The accumulation per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly. By default, the method sums up all attributes per compound, but any function that takes an array and returns an acuumulation over a given axis can be used. For multi-dimensional input arrays, the accumulation is performed along the first axis.

Parameters:
  • attribute (str or array_like) – Attribute or array of values to accumulate. If a numpy.ndarray (or compatible) is provided, its first dimension must have the same length as the total number of atoms in the group.
  • function (callable, optional) – The function performing the accumulation. It must take the array of attribute values to accumulate as its only positional argument and accept an (optional) keyword argument axis allowing to specify the axis along which the accumulation is performed.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the accumulation of all attributes associated with atoms in the group will be returned as a single value. Otherwise, the accumulation of the attributes per Segment, Residue, molecule, or fragment will be returned as a 1d array. Note that, in any case, only the Atoms belonging to the group will be taken into account.
Returns:

Acuumulation of the attribute. If compound is set to 'group', the first dimension of the attribute array will be contracted to a single value. If compound is set to 'segments', 'residues', 'molecules', or 'fragments', the length of the first dimension will correspond to the number of compounds. In all cases, the other dimensions of the returned array will be of the original shape (without the first dimension).

Return type:

float or numpy.ndarray

Raises:
  • ValueError – If the length of a provided attribute array does not correspond to the number of atoms in the group.
  • ValueError – If compound is not one of 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums), or if compound is 'fragments' but the topology doesn’t contain bonds.

Examples

To find the total charge of a given AtomGroup:

>>> sel = u.select_atoms('prop mass > 4.0')
>>> sel.accumulate('charges')

To find the total mass per residue of all CA Atoms:

>>> sel = u.select_atoms('name CA')
>>> sel.accumulate('masses', compound='residues')

To find the maximum atomic charge per fragment of a given AtomGroup:

>>> sel.accumulate('charges', compound="fragments", function=np.max)

New in version 0.20.0.

atoms

An AtomGroup of Atoms present in this ResidueGroup.

The Atoms are ordered locally by Residue in the ResidueGroup. Duplicates are not removed.

bbox(pbc=False)

Return the bounding box of the selection.

The lengths A,B,C of the orthorhombic enclosing box are

L = AtomGroup.bbox()
A,B,C = L[1] - L[0]
Parameters:pbc (bool, optional) – If True, move all Atoms to the primary unit cell before calculation. [False]
Returns:corners – 2x3 array giving corners of bounding box as [[xmin, ymin, zmin], [xmax, ymax, zmax]].
Return type:numpy.ndarray

New in version 0.7.2.

Changed in version 0.8: Added pbc keyword

Changed in version 1.0.0: Removed flags affecting default behaviour

bsphere(pbc=False)

Return the bounding sphere of the selection.

The sphere is calculated relative to the center of geometry.

Parameters:pbc (bool, optional) – If True, move all atoms to the primary unit cell before calculation. [False]
Returns:
  • R (float) – Radius of the bounding sphere.
  • center (numpy.ndarray) – Coordinates of the sphere center as [xcen, ycen, zcen].

New in version 0.7.3.

Changed in version 0.8: Added pbc keyword

center(weights, pbc=False, compound='group', unwrap=False)

Weighted center of (compounds of) the group

Computes the weighted center of Atoms in the group. Weighted centers per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly. If the weights of a compound sum up to zero, the coordinates of that compound’s weighted center will be nan (not a number).

Parameters:
  • weights (array_like or None) – Weights to be used. Setting weights=None is equivalent to passing identical weights for all atoms of the group.
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments', 'residues', 'molecules', or 'fragments', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default [False].
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the weighted center of all atoms in the group will be returned as a single position vector. Else, the weighted centers of each Segment, Residue, molecule, or fragment will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the weighted center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments', 'residues', 'molecules', or 'fragments', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Raises:
  • ValueError – If compound is not one of 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • ValueError – If both ‘pbc’ and ‘unwrap’ set to true.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums) or if compound is 'fragments' but the topology doesn’t contain bonds.

Examples

To find the center of charge of a given AtomGroup:

>>> sel = u.select_atoms('prop mass > 4.0')
>>> sel.center(sel.charges)

To find the centers of mass per residue of all CA Atoms:

>>> sel = u.select_atoms('name CA')
>>> sel.center(sel.masses, compound='residues')

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

center_of_geometry(pbc=False, compound='group', unwrap=False)

Center of geometry of (compounds of) the group.

Computes the center of geometry (a.k.a. centroid) of Atoms in the group. Centers of geometry per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly.

Parameters:
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments' or 'residues', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default False.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the center of geometry of all Atoms in the group will be returned as a single position vector. Else, the centers of geometry of each Segment or Residue will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the geometric center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments' or 'residues', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Changed in version 0.8: Added pbc keyword

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

centroid(pbc=False, compound='group', unwrap=False)

Center of geometry of (compounds of) the group.

Computes the center of geometry (a.k.a. centroid) of Atoms in the group. Centers of geometry per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly.

Parameters:
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments' or 'residues', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default False.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the center of geometry of all Atoms in the group will be returned as a single position vector. Else, the centers of geometry of each Segment or Residue will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the geometric center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments' or 'residues', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Changed in version 0.8: Added pbc keyword

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

concatenate(other)

Concatenate with another Group or Component of the same level.

Duplicate entries and original order is preserved. It is synomymous to the + operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with elements of self and other concatenated
Return type:Group

Example

The order of the original contents (including duplicates) are preserved when performing a concatenation.

>>> ag1 = u.select_atoms('name O')
>>> ag2 = u.select_atoms('name N')
>>> ag3 = ag1 + ag2  # or ag1.concatenate(ag2)
>>> ag3[:3].names
array(['O', 'O', 'O'], dtype=object)
>>> ag3[-3:].names
array(['N', 'N', 'N'], dtype=object)

New in version 0.16.0.

copy()

Get another group identical to this one.

New in version 0.19.0.

difference(other)

Elements from this Group that do not appear in another

This method removes duplicate elements and sorts the result. As such, it is different from subtract(). difference() is synomymous to the - operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements of self that are not in other, without duplicate elements
Return type:Group

New in version 0.16.

dimensions

Obtain a copy of the dimensions of the currently loaded Timestep

groupby(topattrs)

Group together items in this group according to values of topattr

Parameters:topattrs (str or list) – One or more topology attributes to group components by. Single arguments are passed as a string. Multiple arguments are passed as a list.
Returns:Unique values of the multiple combinations of topology attributes as keys, Groups as values.
Return type:dict

Example

To group atoms with the same mass together:

>>> ag.groupby('masses')
{12.010999999999999: <AtomGroup with 462 atoms>,
 14.007: <AtomGroup with 116 atoms>,
 15.999000000000001: <AtomGroup with 134 atoms>}

To group atoms with the same residue name and mass together:

>>> ag.groupby(['resnames', 'masses'])
{('ALA', 1.008): <AtomGroup with 95 atoms>,
 ('ALA', 12.011): <AtomGroup with 57 atoms>,
 ('ALA', 14.007): <AtomGroup with 19 atoms>,
 ('ALA', 15.999): <AtomGroup with 19 atoms>},
 ('ARG', 1.008): <AtomGroup with 169 atoms>,
 ...}
>>> ag.groupby(['resnames', 'masses'])('ALA', 15.999)
 <AtomGroup with 19 atoms>

New in version 0.16.0.

Changed in version 0.18.0: The function accepts multiple attributes

intersection(other)

Group of elements which are in both this Group and another

This method removes duplicate elements and sorts the result. It is synomymous to the & operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the common elements of self and other, without duplicate elements
Return type:Group

Example

Intersections can be used when the select atoms string would become too complicated. For example to find the water atoms which are within 4.0A of two segments:

>>> shell1 = u.select_atoms('resname SOL and around 4.0 segid 1')
>>> shell2 = u.select_atoms('resname SOL and around 4.0 segid 2')
>>> common = shell1 & shell2  # or shell1.intersection(shell2)

See also

union()

New in version 0.16.

is_strict_subset(other)

If this Group is a subset of another Group but not identical

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a strict subset of the other one
Return type:bool

New in version 0.16.

is_strict_superset(other)

If this Group is a superset of another Group but not identical

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a strict superset of the other one
Return type:bool

New in version 0.16.

isdisjoint(other)

If the Group has no elements in common with the other Group

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if the two Groups do not have common elements
Return type:bool

New in version 0.16.

issubset(other)

If all elements of this Group are part of another Group

Note that an empty group is a subset of any group of the same level.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a subset of the other one
Return type:bool

New in version 0.16.

issuperset(other)

If all elements of another Group are part of this Group

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a subset of the other one
Return type:bool

New in version 0.16.

isunique

Boolean indicating whether all components of the group are unique, i.e., the group contains no duplicates.

Examples

>>> ag = u.atoms[[2, 1, 2, 2, 1, 0]]
>>> ag
<AtomGroup with 6 atoms>
>>> ag.isunique
False
>>> ag2 = ag.unique
>>> ag2
<AtomGroup with 3 atoms>
>>> ag2.isunique
True

New in version 0.19.0.

ix

Unique indices of the components in the Group.

ix_array

Unique indices of the components in the Group.

For a Group, ix_array is the same as ix. This method gives a consistent API between components and groups.

See also

ix

n_atoms

Number of Atoms present in this ResidueGroup, including duplicate residues (and thus, duplicate atoms).

Equivalent to len(self.atoms).

n_residues

Number of residues in the ResidueGroup.

Equivalent to len(self).

n_segments

Number of unique segments present in the ResidueGroup.

Equivalent to len(self.segments).

pack_into_box(box=None, inplace=True)

Shift all Atoms in this group to the primary unit cell.

Parameters:
  • box (array_like) – Box dimensions, can be either orthogonal or triclinic information. Cell dimensions must be in an identical to format to those returned by MDAnalysis.coordinates.base.Timestep.dimensions, [lx, ly, lz, alpha, beta, gamma]. If None, uses these timestep dimensions.
  • inplace (bool) – True to change coordinates in place.
Returns:

coords – Shifted atom coordinates.

Return type:

numpy.ndarray

Notes

All atoms will be moved so that they lie between 0 and boxlength \(L_i\) in all dimensions, i.e. the lower left corner of the simulation box is taken to be at (0,0,0):

\[x_i' = x_i - \left\lfloor\frac{x_i}{L_i}\right\rfloor\]

The default is to take unit cell information from the underlying Timestep instance. The optional argument box can be used to provide alternative unit cell information (in the MDAnalysis standard format [Lx, Ly, Lz, alpha, beta, gamma]).

Works with either orthogonal or triclinic box types.

Note

pack_into_box() is identical to wrap() with all default keywords.

New in version 0.8.

residues

The ResidueGroup itself.

See also

copy
return a true copy of the ResidueGroup

Changed in version 0.19.0: In previous versions, this returned a copy, but now the ResidueGroup itself is returned. This should not affect any code but only speed up calculations.

rotate(R, point=(0, 0, 0))

Apply a rotation matrix R to the selection’s coordinates. \(\mathsf{R}\) is a 3x3 orthogonal matrix that transforms a vector \(\mathbf{x} \rightarrow \mathbf{x}'\):

\[\mathbf{x}' = \mathsf{R}\mathbf{x}\]

Atom coordinates are rotated in-place.

Parameters:
  • R (array_like) – 3x3 rotation matrix
  • point (array_like, optional) – Center of rotation
Returns:

Return type:

self

Notes

By default, rotates about the origin point=(0, 0, 0). To rotate a group g around its center of geometry, use g.rotate(R, point=g.center_of_geometry()).

See also

rotateby()
rotate around given axis and angle
MDAnalysis.lib.transformations()
module of all coordinate transforms
rotateby(angle, axis, point=None)

Apply a rotation to the selection’s coordinates.

Parameters:
  • angle (float) – Rotation angle in degrees.
  • axis (array_like) – Rotation axis vector.
  • point (array_like, optional) – Center of rotation. If None then the center of geometry of this group is used.
Returns:

Return type:

self

Notes

The transformation from current coordinates \(\mathbf{x}\) to new coordinates \(\mathbf{x}'\) is

\[\mathbf{x}' = \mathsf{R}\,(\mathbf{x}-\mathbf{p}) + \mathbf{p}\]

where \(\mathsf{R}\) is the rotation by angle around the axis going through point \(\mathbf{p}\).

See also

MDAnalysis.lib.transformations.rotation_matrix()

calculate()
math:mathsf{R}
segments

Get sorted SegmentGroup of the unique segments present in the ResidueGroup.

subtract(other)

Group with elements from this Group that don’t appear in other

The original order of this group is kept, as well as any duplicate elements. If an element of this Group is duplicated and appears in the other Group or Component, then all the occurences of that element are removed from the returned Group.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements of self that are not in other, conserves order and duplicates.
Return type:Group

Example

Unlike difference() this method will not sort or remove duplicates.

>>> ag1 = u.atoms[[3, 3, 2, 2, 1, 1]]
>>> ag2 = u.atoms[2]
>>> ag3 = ag1 - ag2  # or ag1.subtract(ag2)
>>> ag1.indices
array([3, 3, 1, 1])

New in version 0.16.

symmetric_difference(other)

Group of elements which are only in one of this Group or another

This method removes duplicate elements and the result is sorted. It is synomym to the ^ operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements that are in self or in other but not in both, without duplicate elements
Return type:Group

Example

>>> ag1 = u.atoms[[0, 1, 5, 3, 3, 2]]
>>> ag2 = u.atoms[[4, 4, 6, 2, 3, 5]]
>>> ag3 = ag1 ^ ag2  # or ag1.symmetric_difference(ag2)
>>> ag3.indices  # 0 and 1 are only in ag1, 4 and 6 are only in ag2
[0, 1, 4, 6]

See also

difference()

New in version 0.16.

transform(M)

Apply homogenous transformation matrix M to the coordinates.

Atom coordinates are rotated and translated in-place.

Parameters:M (array_like) – 4x4 matrix with the rotation in R = M[:3, :3] and the translation in t = M[:3, 3].
Returns:
Return type:self

See also

MDAnalysis.lib.transformations()
module of all coordinate transforms

Notes

The rotation \(\mathsf{R}\) is about the origin and is applied before the translation \(\mathbf{t}\):

\[\mathbf{x}' = \mathsf{R}\mathbf{x} + \mathbf{t}\]
translate(t)

Apply translation vector t to the selection’s coordinates.

Atom coordinates are translated in-place.

Parameters:t (array_like) – vector to translate coordinates with
Returns:
Return type:self

See also

MDAnalysis.lib.transformations()
module of all coordinate transforms

Notes

The method applies a translation to the AtomGroup from current coordinates \(\mathbf{x}\) to new coordinates \(\mathbf{x}'\):

\[\mathbf{x}' = \mathbf{x} + \mathbf{t}\]
union(other)

Group of elements either in this Group or another

On the contrary to concatenation, this method sort the elements and removes duplicate ones. It is synomymous to the | operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the combined elements of self and other, without duplicate elements
Return type:Group

Example

In contrast to concatenate(), any duplicates are dropped and the result is sorted.

>>> ag1 = u.select_atoms('name O')
>>> ag2 = u.select_atoms('name N')
>>> ag3 = ag1 | ag2  # or ag1.union(ag2)
>>> ag3[:3].names
array(['N', 'O', 'N'], dtype=object)

New in version 0.16.

unique

Return a ResidueGroup containing sorted and unique Residues only.

If the ResidueGroup is unique, this is the group itself.

Examples

>>> rg = u.residues[[2, 1, 2, 2, 1, 0]]
>>> rg
<ResidueGroup with 6 residues>
>>> rg.ix
array([2, 1, 2, 2, 1, 0])
>>> rg2 = rg.unique
>>> rg2
<ResidueGroup with 3 residues>
>>> rg2.ix
array([0, 1, 2])
>>> rg2.unique is rg2
True

New in version 0.16.0.

Changed in version 0.19.0: If the ResidueGroup is already unique, ResidueGroup.unique now returns the group itself instead of a copy.

universe

The underlying Universe the group belongs to.

unwrap(compound='fragments', reference='com', inplace=True)

Move atoms of this group so that bonds within the group’s compounds aren’t split across periodic boundaries.

This function is most useful when atoms have been packed into the primary unit cell, causing breaks mid-molecule, with the molecule then appearing on either side of the unit cell. This is problematic for operations such as calculating the center of mass of the molecule.

+-----------+       +-----------+
|           |       |           |
| 6       3 |       |         3 | 6
| !       ! |       |         ! | !
|-5-8   1-2-|  ==>  |       1-2-|-5-8
| !       ! |       |         ! | !
| 7       4 |       |         4 | 7
|           |       |           |
+-----------+       +-----------+
Parameters:
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – Which type of component to unwrap. Note that, in any case, all atoms within each compound must be interconnected by bonds, i.e., compounds must correspond to (parts of) molecules.
  • reference ({'com', 'cog', None}, optional) – If 'com' (center of mass) or 'cog' (center of geometry), the unwrapped compounds will be shifted so that their individual reference point lies within the primary unit cell. If None, no such shift is performed.
  • inplace (bool, optional) – If True, coordinates are modified in place.
Returns:

coords – Unwrapped atom coordinate array of shape (n, 3).

Return type:

numpy.ndarray

Raises:
  • NoDataError – If compound is 'molecules' but the underlying topology does not contain molecule information, or if reference is 'com' but the topology does not contain masses.
  • ValueError – If reference is not one of 'com', 'cog', or None, or if reference is 'com' and the total mass of any compound is zero.

Note

Be aware of the fact that only atoms belonging to the group will be unwrapped! If you want entire molecules to be unwrapped, make sure that all atoms of these molecules are part of the group.

An AtomGroup containing all atoms of all fragments in the group ag can be created with:

all_frag_atoms = sum(ag.fragments)

See also

make_whole(), wrap(), pack_into_box(), apply_PBC()

New in version 0.20.0.

wrap(compound='atoms', center='com', box=None, inplace=True)

Shift the contents of this group back into the primary unit cell according to periodic boundary conditions.

Specifying a compound will keep the Atoms in each compound together during the process. If compound is different from 'atoms', each compound as a whole will be shifted so that its center lies within the primary unit cell.

Parameters:
  • compound ({'atoms', 'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – Which type of component to keep together during wrapping. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • center ({'com', 'cog'}) – How to define the center of a given group of atoms. If compound is 'atoms', this parameter is meaningless and therefore ignored.
  • box (array_like, optional) –

    The unitcell dimensions of the system, which can be orthogonal or triclinic and must be provided in the same format as returned by MDAnalysis.coordinates.base.Timestep.dimensions:

    [lx, ly, lz, alpha, beta, gamma].

    If None, uses the dimensions of the current time step.

  • inplace (bool, optional) – If True, coordinates will be changed in place.
Returns:

Array of wrapped atom coordinates of dtype np.float32 and shape (len(self.atoms.n_atoms), 3)

Return type:

numpy.ndarray

Raises:
  • ValueError – If compound is not one of 'atoms', 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums) or if compound is 'fragments' but the topology doesn’t contain bonds or if center is 'com' but the topology doesn’t contain masses.

Notes

All atoms of the group will be moved so that the centers of its compounds lie within the primary periodic image. For orthorhombic unit cells, the primary periodic image is defined as the half-open interval \([0,L_i)\) between \(0\) and boxlength \(L_i\) in all dimensions \(i\in\{x,y,z\}\), i.e., the origin of the of the simulation box is taken to be at the origin \((0,0,0)\) of the euclidian coordinate system. A compound center residing at position \(x_i\) in dimension \(i\) will be shifted to \(x_i'\) according to

\[x_i' = x_i - \left\lfloor\frac{x_i}{L_i}\right\rfloor\,.\]

When specifying a compound, the translation is calculated based on each compound. The same translation is applied to all atoms within this compound, meaning it will not be broken by the shift. This might however mean that not all atoms of a compound will be inside the unit cell after wrapping, but rather will be the center of the compound.

Be aware of the fact that only atoms belonging to the group will be taken into account!

center allows to define how the center of each group is computed. This can be either 'com' for center of mass, or 'cog' for center of geometry.

box allows a unit cell to be given for the transformation. If not specified, the dimensions information from the current Timestep will be used.

Note

AtomGroup.wrap() is currently faster than ResidueGroup.wrap() or SegmentGroup.wrap().

See also

pack_into_box(), unwrap(), MDanalysis.lib.distances.apply_PBC()

New in version 0.9.2.

Changed in version 0.20.0: The method only acts on atoms belonging to the group and returns the wrapped positions as a numpy.ndarray. Added optional argument inplace.

class MDAnalysis.core.groups.SegmentGroup(*args)[source]

SegmentGroup base class.

This class is used by a Universe for generating its Topology-specific SegmentGroup class. All the TopologyAttr components are obtained from GroupBase, so this class only includes ad-hoc methods specific to SegmentGroups.

SegmentGroups can be compared and combined using group operators. See the list of these operators on GroupBase.

Deprecated since version 0.16.2: Instant selectors of Segments will be removed in the 1.0 release.

Changed in version 1.0.0: Removed instant selectors, use select_atoms instead

accumulate(attribute, function=<function sum>, compound='group')

Accumulates the attribute associated with (compounds of) the group.

Accumulates the attribute of Atoms in the group. The accumulation per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly. By default, the method sums up all attributes per compound, but any function that takes an array and returns an acuumulation over a given axis can be used. For multi-dimensional input arrays, the accumulation is performed along the first axis.

Parameters:
  • attribute (str or array_like) – Attribute or array of values to accumulate. If a numpy.ndarray (or compatible) is provided, its first dimension must have the same length as the total number of atoms in the group.
  • function (callable, optional) – The function performing the accumulation. It must take the array of attribute values to accumulate as its only positional argument and accept an (optional) keyword argument axis allowing to specify the axis along which the accumulation is performed.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the accumulation of all attributes associated with atoms in the group will be returned as a single value. Otherwise, the accumulation of the attributes per Segment, Residue, molecule, or fragment will be returned as a 1d array. Note that, in any case, only the Atoms belonging to the group will be taken into account.
Returns:

Acuumulation of the attribute. If compound is set to 'group', the first dimension of the attribute array will be contracted to a single value. If compound is set to 'segments', 'residues', 'molecules', or 'fragments', the length of the first dimension will correspond to the number of compounds. In all cases, the other dimensions of the returned array will be of the original shape (without the first dimension).

Return type:

float or numpy.ndarray

Raises:
  • ValueError – If the length of a provided attribute array does not correspond to the number of atoms in the group.
  • ValueError – If compound is not one of 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums), or if compound is 'fragments' but the topology doesn’t contain bonds.

Examples

To find the total charge of a given AtomGroup:

>>> sel = u.select_atoms('prop mass > 4.0')
>>> sel.accumulate('charges')

To find the total mass per residue of all CA Atoms:

>>> sel = u.select_atoms('name CA')
>>> sel.accumulate('masses', compound='residues')

To find the maximum atomic charge per fragment of a given AtomGroup:

>>> sel.accumulate('charges', compound="fragments", function=np.max)

New in version 0.20.0.

atoms

An AtomGroup of Atoms present in this SegmentGroup.

The Atoms are ordered locally by Residue, which are further ordered by Segment in the SegmentGroup. Duplicates are not removed.

bbox(pbc=False)

Return the bounding box of the selection.

The lengths A,B,C of the orthorhombic enclosing box are

L = AtomGroup.bbox()
A,B,C = L[1] - L[0]
Parameters:pbc (bool, optional) – If True, move all Atoms to the primary unit cell before calculation. [False]
Returns:corners – 2x3 array giving corners of bounding box as [[xmin, ymin, zmin], [xmax, ymax, zmax]].
Return type:numpy.ndarray

New in version 0.7.2.

Changed in version 0.8: Added pbc keyword

Changed in version 1.0.0: Removed flags affecting default behaviour

bsphere(pbc=False)

Return the bounding sphere of the selection.

The sphere is calculated relative to the center of geometry.

Parameters:pbc (bool, optional) – If True, move all atoms to the primary unit cell before calculation. [False]
Returns:
  • R (float) – Radius of the bounding sphere.
  • center (numpy.ndarray) – Coordinates of the sphere center as [xcen, ycen, zcen].

New in version 0.7.3.

Changed in version 0.8: Added pbc keyword

center(weights, pbc=False, compound='group', unwrap=False)

Weighted center of (compounds of) the group

Computes the weighted center of Atoms in the group. Weighted centers per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly. If the weights of a compound sum up to zero, the coordinates of that compound’s weighted center will be nan (not a number).

Parameters:
  • weights (array_like or None) – Weights to be used. Setting weights=None is equivalent to passing identical weights for all atoms of the group.
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments', 'residues', 'molecules', or 'fragments', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default [False].
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the weighted center of all atoms in the group will be returned as a single position vector. Else, the weighted centers of each Segment, Residue, molecule, or fragment will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the weighted center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments', 'residues', 'molecules', or 'fragments', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Raises:
  • ValueError – If compound is not one of 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • ValueError – If both ‘pbc’ and ‘unwrap’ set to true.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums) or if compound is 'fragments' but the topology doesn’t contain bonds.

Examples

To find the center of charge of a given AtomGroup:

>>> sel = u.select_atoms('prop mass > 4.0')
>>> sel.center(sel.charges)

To find the centers of mass per residue of all CA Atoms:

>>> sel = u.select_atoms('name CA')
>>> sel.center(sel.masses, compound='residues')

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

center_of_geometry(pbc=False, compound='group', unwrap=False)

Center of geometry of (compounds of) the group.

Computes the center of geometry (a.k.a. centroid) of Atoms in the group. Centers of geometry per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly.

Parameters:
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments' or 'residues', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default False.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the center of geometry of all Atoms in the group will be returned as a single position vector. Else, the centers of geometry of each Segment or Residue will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the geometric center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments' or 'residues', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Changed in version 0.8: Added pbc keyword

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

centroid(pbc=False, compound='group', unwrap=False)

Center of geometry of (compounds of) the group.

Computes the center of geometry (a.k.a. centroid) of Atoms in the group. Centers of geometry per Residue, Segment, molecule, or fragment can be obtained by setting the compound parameter accordingly.

Parameters:
  • pbc (bool, optional) – If True and compound is 'group', move all atoms to the primary unit cell before calculation. If True and compound is 'segments' or 'residues', the center of each compound will be calculated without moving any Atoms to keep the compounds intact. Instead, the resulting position vectors will be moved to the primary unit cell after calculation. Default False.
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – If 'group', the center of geometry of all Atoms in the group will be returned as a single position vector. Else, the centers of geometry of each Segment or Residue will be returned as an array of position vectors, i.e. a 2d array. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • unwrap (bool, optional) – If True, compounds will be unwrapped before computing their centers.
Returns:

center – Position vector(s) of the geometric center(s) of the group. If compound was set to 'group', the output will be a single position vector. If compound was set to 'segments' or 'residues', the output will be a 2d array of shape (n, 3) where n is the number of compounds.

Return type:

numpy.ndarray

Changed in version 0.8: Added pbc keyword

Changed in version 0.19.0: Added compound parameter

Changed in version 0.20.0: Added 'molecules' and 'fragments' compounds

Changed in version 0.20.0: Added unwrap parameter

Changed in version 1.0.0: Removed flags affecting default behaviour

concatenate(other)

Concatenate with another Group or Component of the same level.

Duplicate entries and original order is preserved. It is synomymous to the + operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with elements of self and other concatenated
Return type:Group

Example

The order of the original contents (including duplicates) are preserved when performing a concatenation.

>>> ag1 = u.select_atoms('name O')
>>> ag2 = u.select_atoms('name N')
>>> ag3 = ag1 + ag2  # or ag1.concatenate(ag2)
>>> ag3[:3].names
array(['O', 'O', 'O'], dtype=object)
>>> ag3[-3:].names
array(['N', 'N', 'N'], dtype=object)

New in version 0.16.0.

copy()

Get another group identical to this one.

New in version 0.19.0.

difference(other)

Elements from this Group that do not appear in another

This method removes duplicate elements and sorts the result. As such, it is different from subtract(). difference() is synomymous to the - operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements of self that are not in other, without duplicate elements
Return type:Group

New in version 0.16.

dimensions

Obtain a copy of the dimensions of the currently loaded Timestep

groupby(topattrs)

Group together items in this group according to values of topattr

Parameters:topattrs (str or list) – One or more topology attributes to group components by. Single arguments are passed as a string. Multiple arguments are passed as a list.
Returns:Unique values of the multiple combinations of topology attributes as keys, Groups as values.
Return type:dict

Example

To group atoms with the same mass together:

>>> ag.groupby('masses')
{12.010999999999999: <AtomGroup with 462 atoms>,
 14.007: <AtomGroup with 116 atoms>,
 15.999000000000001: <AtomGroup with 134 atoms>}

To group atoms with the same residue name and mass together:

>>> ag.groupby(['resnames', 'masses'])
{('ALA', 1.008): <AtomGroup with 95 atoms>,
 ('ALA', 12.011): <AtomGroup with 57 atoms>,
 ('ALA', 14.007): <AtomGroup with 19 atoms>,
 ('ALA', 15.999): <AtomGroup with 19 atoms>},
 ('ARG', 1.008): <AtomGroup with 169 atoms>,
 ...}
>>> ag.groupby(['resnames', 'masses'])('ALA', 15.999)
 <AtomGroup with 19 atoms>

New in version 0.16.0.

Changed in version 0.18.0: The function accepts multiple attributes

intersection(other)

Group of elements which are in both this Group and another

This method removes duplicate elements and sorts the result. It is synomymous to the & operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the common elements of self and other, without duplicate elements
Return type:Group

Example

Intersections can be used when the select atoms string would become too complicated. For example to find the water atoms which are within 4.0A of two segments:

>>> shell1 = u.select_atoms('resname SOL and around 4.0 segid 1')
>>> shell2 = u.select_atoms('resname SOL and around 4.0 segid 2')
>>> common = shell1 & shell2  # or shell1.intersection(shell2)

See also

union()

New in version 0.16.

is_strict_subset(other)

If this Group is a subset of another Group but not identical

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a strict subset of the other one
Return type:bool

New in version 0.16.

is_strict_superset(other)

If this Group is a superset of another Group but not identical

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a strict superset of the other one
Return type:bool

New in version 0.16.

isdisjoint(other)

If the Group has no elements in common with the other Group

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if the two Groups do not have common elements
Return type:bool

New in version 0.16.

issubset(other)

If all elements of this Group are part of another Group

Note that an empty group is a subset of any group of the same level.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a subset of the other one
Return type:bool

New in version 0.16.

issuperset(other)

If all elements of another Group are part of this Group

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:True if this Group is a subset of the other one
Return type:bool

New in version 0.16.

isunique

Boolean indicating whether all components of the group are unique, i.e., the group contains no duplicates.

Examples

>>> ag = u.atoms[[2, 1, 2, 2, 1, 0]]
>>> ag
<AtomGroup with 6 atoms>
>>> ag.isunique
False
>>> ag2 = ag.unique
>>> ag2
<AtomGroup with 3 atoms>
>>> ag2.isunique
True

New in version 0.19.0.

ix

Unique indices of the components in the Group.

ix_array

Unique indices of the components in the Group.

For a Group, ix_array is the same as ix. This method gives a consistent API between components and groups.

See also

ix

n_atoms

Number of atoms present in the SegmentGroup, including duplicate segments (and thus, duplicate atoms).

Equivalent to len(self.atoms).

n_residues

Number of residues present in this SegmentGroup, including duplicate segments (and thus, residues).

Equivalent to len(self.residues).

n_segments

Number of segments in the SegmentGroup.

Equivalent to len(self).

pack_into_box(box=None, inplace=True)

Shift all Atoms in this group to the primary unit cell.

Parameters:
  • box (array_like) – Box dimensions, can be either orthogonal or triclinic information. Cell dimensions must be in an identical to format to those returned by MDAnalysis.coordinates.base.Timestep.dimensions, [lx, ly, lz, alpha, beta, gamma]. If None, uses these timestep dimensions.
  • inplace (bool) – True to change coordinates in place.
Returns:

coords – Shifted atom coordinates.

Return type:

numpy.ndarray

Notes

All atoms will be moved so that they lie between 0 and boxlength \(L_i\) in all dimensions, i.e. the lower left corner of the simulation box is taken to be at (0,0,0):

\[x_i' = x_i - \left\lfloor\frac{x_i}{L_i}\right\rfloor\]

The default is to take unit cell information from the underlying Timestep instance. The optional argument box can be used to provide alternative unit cell information (in the MDAnalysis standard format [Lx, Ly, Lz, alpha, beta, gamma]).

Works with either orthogonal or triclinic box types.

Note

pack_into_box() is identical to wrap() with all default keywords.

New in version 0.8.

residues

A ResidueGroup of Residues present in this SegmentGroup.

The Residues are ordered locally by Segment in the SegmentGroup. Duplicates are not removed.

rotate(R, point=(0, 0, 0))

Apply a rotation matrix R to the selection’s coordinates. \(\mathsf{R}\) is a 3x3 orthogonal matrix that transforms a vector \(\mathbf{x} \rightarrow \mathbf{x}'\):

\[\mathbf{x}' = \mathsf{R}\mathbf{x}\]

Atom coordinates are rotated in-place.

Parameters:
  • R (array_like) – 3x3 rotation matrix
  • point (array_like, optional) – Center of rotation
Returns:

Return type:

self

Notes

By default, rotates about the origin point=(0, 0, 0). To rotate a group g around its center of geometry, use g.rotate(R, point=g.center_of_geometry()).

See also

rotateby()
rotate around given axis and angle
MDAnalysis.lib.transformations()
module of all coordinate transforms
rotateby(angle, axis, point=None)

Apply a rotation to the selection’s coordinates.

Parameters:
  • angle (float) – Rotation angle in degrees.
  • axis (array_like) – Rotation axis vector.
  • point (array_like, optional) – Center of rotation. If None then the center of geometry of this group is used.
Returns:

Return type:

self

Notes

The transformation from current coordinates \(\mathbf{x}\) to new coordinates \(\mathbf{x}'\) is

\[\mathbf{x}' = \mathsf{R}\,(\mathbf{x}-\mathbf{p}) + \mathbf{p}\]

where \(\mathsf{R}\) is the rotation by angle around the axis going through point \(\mathbf{p}\).

See also

MDAnalysis.lib.transformations.rotation_matrix()

calculate()
math:mathsf{R}
segments

The SegmentGroup itself.

See also

copy
return a true copy of the SegmentGroup

Changed in version 0.19.0: In previous versions, this returned a copy, but now the SegmentGroup itself is returned. This should not affect any code but only speed up calculations.

subtract(other)

Group with elements from this Group that don’t appear in other

The original order of this group is kept, as well as any duplicate elements. If an element of this Group is duplicated and appears in the other Group or Component, then all the occurences of that element are removed from the returned Group.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements of self that are not in other, conserves order and duplicates.
Return type:Group

Example

Unlike difference() this method will not sort or remove duplicates.

>>> ag1 = u.atoms[[3, 3, 2, 2, 1, 1]]
>>> ag2 = u.atoms[2]
>>> ag3 = ag1 - ag2  # or ag1.subtract(ag2)
>>> ag1.indices
array([3, 3, 1, 1])

New in version 0.16.

symmetric_difference(other)

Group of elements which are only in one of this Group or another

This method removes duplicate elements and the result is sorted. It is synomym to the ^ operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the elements that are in self or in other but not in both, without duplicate elements
Return type:Group

Example

>>> ag1 = u.atoms[[0, 1, 5, 3, 3, 2]]
>>> ag2 = u.atoms[[4, 4, 6, 2, 3, 5]]
>>> ag3 = ag1 ^ ag2  # or ag1.symmetric_difference(ag2)
>>> ag3.indices  # 0 and 1 are only in ag1, 4 and 6 are only in ag2
[0, 1, 4, 6]

See also

difference()

New in version 0.16.

transform(M)

Apply homogenous transformation matrix M to the coordinates.

Atom coordinates are rotated and translated in-place.

Parameters:M (array_like) – 4x4 matrix with the rotation in R = M[:3, :3] and the translation in t = M[:3, 3].
Returns:
Return type:self

See also

MDAnalysis.lib.transformations()
module of all coordinate transforms

Notes

The rotation \(\mathsf{R}\) is about the origin and is applied before the translation \(\mathbf{t}\):

\[\mathbf{x}' = \mathsf{R}\mathbf{x} + \mathbf{t}\]
translate(t)

Apply translation vector t to the selection’s coordinates.

Atom coordinates are translated in-place.

Parameters:t (array_like) – vector to translate coordinates with
Returns:
Return type:self

See also

MDAnalysis.lib.transformations()
module of all coordinate transforms

Notes

The method applies a translation to the AtomGroup from current coordinates \(\mathbf{x}\) to new coordinates \(\mathbf{x}'\):

\[\mathbf{x}' = \mathbf{x} + \mathbf{t}\]
union(other)

Group of elements either in this Group or another

On the contrary to concatenation, this method sort the elements and removes duplicate ones. It is synomymous to the | operator.

Parameters:other (Group or Component) – Group or Component with other.level same as self.level
Returns:Group with the combined elements of self and other, without duplicate elements
Return type:Group

Example

In contrast to concatenate(), any duplicates are dropped and the result is sorted.

>>> ag1 = u.select_atoms('name O')
>>> ag2 = u.select_atoms('name N')
>>> ag3 = ag1 | ag2  # or ag1.union(ag2)
>>> ag3[:3].names
array(['N', 'O', 'N'], dtype=object)

New in version 0.16.

unique

Return a SegmentGroup containing sorted and unique Segments only.

If the SegmentGroup is unique, this is the group itself.

Examples

>>> sg = u.segments[[2, 1, 2, 2, 1, 0]]
>>> sg
<SegmentGroup with 6 segments>
>>> sg.ix
array([2, 1, 2, 2, 1, 0])
>>> sg2 = sg.unique
>>> sg2
<SegmentGroup with 3 segments>
>>> sg2.ix
array([0, 1, 2])
>>> sg2.unique is sg2
True

New in version 0.16.0.

Changed in version 0.19.0: If the SegmentGroup is already unique, SegmentGroup.unique now returns the group itself instead of a copy.

universe

The underlying Universe the group belongs to.

unwrap(compound='fragments', reference='com', inplace=True)

Move atoms of this group so that bonds within the group’s compounds aren’t split across periodic boundaries.

This function is most useful when atoms have been packed into the primary unit cell, causing breaks mid-molecule, with the molecule then appearing on either side of the unit cell. This is problematic for operations such as calculating the center of mass of the molecule.

+-----------+       +-----------+
|           |       |           |
| 6       3 |       |         3 | 6
| !       ! |       |         ! | !
|-5-8   1-2-|  ==>  |       1-2-|-5-8
| !       ! |       |         ! | !
| 7       4 |       |         4 | 7
|           |       |           |
+-----------+       +-----------+
Parameters:
  • compound ({'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – Which type of component to unwrap. Note that, in any case, all atoms within each compound must be interconnected by bonds, i.e., compounds must correspond to (parts of) molecules.
  • reference ({'com', 'cog', None}, optional) – If 'com' (center of mass) or 'cog' (center of geometry), the unwrapped compounds will be shifted so that their individual reference point lies within the primary unit cell. If None, no such shift is performed.
  • inplace (bool, optional) – If True, coordinates are modified in place.
Returns:

coords – Unwrapped atom coordinate array of shape (n, 3).

Return type:

numpy.ndarray

Raises:
  • NoDataError – If compound is 'molecules' but the underlying topology does not contain molecule information, or if reference is 'com' but the topology does not contain masses.
  • ValueError – If reference is not one of 'com', 'cog', or None, or if reference is 'com' and the total mass of any compound is zero.

Note

Be aware of the fact that only atoms belonging to the group will be unwrapped! If you want entire molecules to be unwrapped, make sure that all atoms of these molecules are part of the group.

An AtomGroup containing all atoms of all fragments in the group ag can be created with:

all_frag_atoms = sum(ag.fragments)

See also

make_whole(), wrap(), pack_into_box(), apply_PBC()

New in version 0.20.0.

wrap(compound='atoms', center='com', box=None, inplace=True)

Shift the contents of this group back into the primary unit cell according to periodic boundary conditions.

Specifying a compound will keep the Atoms in each compound together during the process. If compound is different from 'atoms', each compound as a whole will be shifted so that its center lies within the primary unit cell.

Parameters:
  • compound ({'atoms', 'group', 'segments', 'residues', 'molecules', 'fragments'}, optional) – Which type of component to keep together during wrapping. Note that, in any case, only the positions of Atoms belonging to the group will be taken into account.
  • center ({'com', 'cog'}) – How to define the center of a given group of atoms. If compound is 'atoms', this parameter is meaningless and therefore ignored.
  • box (array_like, optional) –

    The unitcell dimensions of the system, which can be orthogonal or triclinic and must be provided in the same format as returned by MDAnalysis.coordinates.base.Timestep.dimensions:

    [lx, ly, lz, alpha, beta, gamma].

    If None, uses the dimensions of the current time step.

  • inplace (bool, optional) – If True, coordinates will be changed in place.
Returns:

Array of wrapped atom coordinates of dtype np.float32 and shape (len(self.atoms.n_atoms), 3)

Return type:

numpy.ndarray

Raises:
  • ValueError – If compound is not one of 'atoms', 'group', 'segments', 'residues', 'molecules', or 'fragments'.
  • NoDataError – If compound is 'molecule' but the topology doesn’t contain molecule information (molnums) or if compound is 'fragments' but the topology doesn’t contain bonds or if center is 'com' but the topology doesn’t contain masses.

Notes

All atoms of the group will be moved so that the centers of its compounds lie within the primary periodic image. For orthorhombic unit cells, the primary periodic image is defined as the half-open interval \([0,L_i)\) between \(0\) and boxlength \(L_i\) in all dimensions \(i\in\{x,y,z\}\), i.e., the origin of the of the simulation box is taken to be at the origin \((0,0,0)\) of the euclidian coordinate system. A compound center residing at position \(x_i\) in dimension \(i\) will be shifted to \(x_i'\) according to

\[x_i' = x_i - \left\lfloor\frac{x_i}{L_i}\right\rfloor\,.\]

When specifying a compound, the translation is calculated based on each compound. The same translation is applied to all atoms within this compound, meaning it will not be broken by the shift. This might however mean that not all atoms of a compound will be inside the unit cell after wrapping, but rather will be the center of the compound.

Be aware of the fact that only atoms belonging to the group will be taken into account!

center allows to define how the center of each group is computed. This can be either 'com' for center of mass, or 'cog' for center of geometry.

box allows a unit cell to be given for the transformation. If not specified, the dimensions information from the current Timestep will be used.

Note

AtomGroup.wrap() is currently faster than ResidueGroup.wrap() or SegmentGroup.wrap().

See also

pack_into_box(), unwrap(), MDanalysis.lib.distances.apply_PBC()

New in version 0.9.2.

Changed in version 0.20.0: The method only acts on atoms belonging to the group and returns the wrapped positions as a numpy.ndarray. Added optional argument inplace.

class MDAnalysis.core.groups.UpdatingAtomGroup(base_group, selections, strings)[source]

AtomGroup subclass that dynamically updates its selected atoms.

Accessing any attribute/method of an UpdatingAtomGroup instance triggers a check for the last frame the group was updated. If the last frame matches the current trajectory frame, the attribute is returned normally; otherwise the group is updated (the stored selections are re-applied), and only then is the attribute returned.

New in version 0.16.0.

Parameters:
  • base_group (AtomGroup) – group on which selections are to be applied.
  • selections (a tuple of Selection) – instances selections ready to be applied to base_group.
atoms

Get a static AtomGroup identical to the group of currently selected Atoms in the UpdatingAtomGroup.

By returning a static AtomGroup it becomes possible to compare the contents of the group between trajectory frames. See the Example below.

Note

The atoms attribute of an UpdatingAtomGroup behaves differently from AtomGroup.atoms: the latter returns the AtomGroup itself whereas the former returns a AtomGroup and not an UpdatingAtomGroup (for this, use UpdatingAtomGroup.copy()).

Example

The static atoms allows comparison of groups of atoms between frames. For example, track water molecules that move in and out of a solvation shell of a protein:

u = mda.Universe(TPR, XTC)
water_shell = u.select_atoms("name OW and around 3.5 protein", updating=True)
water_shell_prev = water_shell.atoms

for ts in u.trajectory:
    exchanged = water_shell - water_shell_prev

    print(ts.time, "waters in shell =", water_shell.n_residues)
    print(ts.time, "waters that exchanged = ", exchanged.n_residues)
    print(ts.time, "waters that remained bound = ",
          water_shell.n_residues - exchanged.n_residues)

    water_shell_prev = water_shell.atoms

By remembering the atoms of the current time step in water_shell_prev, it becomes possible to use the subtraction of AtomGroups to find the water molecules that changed.

See also

copy
return a true copy of the UpdatingAtomGroup
copy()[source]

Get another UpdatingAtomGroup identical to this one.

New in version 0.19.0.

is_uptodate

Checks whether the selection needs updating based on frame number only.

Modifications to the coordinate data that render selections stale are not caught, and in those cases is_uptodate may return an erroneous value.

Returns:True if the group’s selection is up-to-date, False otherwise.
Return type:bool
update_selection()[source]

Forces the reevaluation and application of the group’s selection(s).

This method is triggered automatically when accessing attributes, if the last update occurred under a different trajectory frame.

11.1.2.1.2. Chemical units

class MDAnalysis.core.groups.Atom(ix, u)[source]

Atom base class.

This class is used by a Universe for generating its Topology-specific Atom class. All the TopologyAttr components are obtained from ComponentBase, so this class only includes ad-hoc methods specific to Atoms.

force

Force on the atom.

The force can be changed by assigning an array of shape (3,).

Note

changing the force is not reflected in any files; reading any frame from the trajectory will replace the change with that from the file

Raises:NoDataError – If the underlying Timestep does not contain forces.
ix

Unique index of this component.

If this component is an Atom, this is the index of the Atom. If it is a Residue, this is the index of the Residue. If it is a Segment, this is the index of the Segment.

ix_array

Unique index of this component as an array.

This method gives a consistent API between components and groups.

See also

ix

position

Coordinates of the atom.

The position can be changed by assigning an array of length (3,).

Note

changing the position is not reflected in any files; reading any frame from the trajectory will replace the change with that from the file

Raises:NoDataError – If the underlying Timestep does not contain positions.
velocity

Velocity of the atom.

The velocity can be changed by assigning an array of shape (3,).

Note

changing the velocity is not reflected in any files; reading any frame from the trajectory will replace the change with that from the file

Raises:NoDataError – If the underlying Timestep does not contain velocities.
class MDAnalysis.core.groups.Residue(ix, u)[source]

Residue base class.

This class is used by a Universe for generating its Topology-specific Residue class. All the TopologyAttr components are obtained from ComponentBase, so this class only includes ad-hoc methods specific to Residues.

atoms

An AtomGroup of Atoms present in this Residue.

ix

Unique index of this component.

If this component is an Atom, this is the index of the Atom. If it is a Residue, this is the index of the Residue. If it is a Segment, this is the index of the Segment.

ix_array

Unique index of this component as an array.

This method gives a consistent API between components and groups.

See also

ix

segment

The Segment this Residue belongs to.

class MDAnalysis.core.groups.Segment(ix, u)[source]

Segment base class.

This class is used by a Universe for generating its Topology-specific Segment class. All the TopologyAttr components are obtained from ComponentBase, so this class only includes ad-hoc methods specific to Segments.

Deprecated since version 0.16.2: Instant selectors of Segments will be removed in the 1.0 release.

Changed in version 1.0.0: Removed instant selectors, use either segment.residues[…] to select residue by number, or segment.residues[segment.residue.resnames = …] to select by resname.

atoms

An AtomGroup of Atoms present in this Segment.

ix

Unique index of this component.

If this component is an Atom, this is the index of the Atom. If it is a Residue, this is the index of the Residue. If it is a Segment, this is the index of the Segment.

ix_array

Unique index of this component as an array.

This method gives a consistent API between components and groups.

See also

ix

residues

A ResidueGroup of Residues present in this Segment.

11.1.2.1.3. Levels

Each of the above classes has a level attribute. This can be used to verify that two objects are of the same level, or to access a particular class:

u = mda.Universe()

ag = u.atoms[:10]
at = u.atoms[11]

ag.level == at.level  # Returns True

ag.level.singular  # Returns Atom class
at.level.plural  # Returns AtomGroup class