7.2. Fundamental building blocks — MDAnalysis.core.AtomGroup

The most important data structure in MDAnalysis is the AtomGroup, which contains Atom instances.

A Universe is the user-visible entry point and collects all information needed to analyze a structure or a whole trajectory.

Segments and residues are a way to refer to a collection of atoms. By convention, a Residue is a single amino acid, or a water molecule, ion, or ligand. A Segment is a collection of residues such as a whole protein or a chain in a protein or all the water in the system.

7.2.1. Class Hierarchy

A Universe contains Segments, which contain Residues, which contain Atoms; all containers are derived from AtomGroup, and thus one can always analyze them as a collection of atoms, independent of the hierarchical level.

Each Atom can only belong to a single Residue, and a Residue belongs to one specific Segment. This hierarchy can be described as

Segment > Residue > Atom

Depending on the use case, it can be more convenient to access data on, for instance, the basis of residues than atoms, or to write out individual chains (segments) of a protein. MDAnalysis simply provides access to these groupings and keeps track of where an atom belongs. Each object provides three attributes (atoms, residues or residue, segments or segment) that give access to the tiers in the hierarchy that the object belongs to.

7.2.2. Manipulating atoms, residues, and segments

When working with MDAnalysis it is useful to remember that the fundamental object is the Atom. Each particle in the topology is represented by exactly one Atom instance. One Atom, however, can be a member of multiple AtomGroup collections, for instance from different selections even though they all refer to the same Atom object. Thus, changing a property of a specific and Atom in one AtomGroup changes it “everywhere”.

The same is mostly true for Residue instances although they are derived from Atom instances: all Atom objects with the same Atom.resid are bundled into a single Residue with Residue.id = resid. This means that just changing, say, the residue name with a command such as

>>> r = u.select_atoms("resid 99").residues[0]
>>> print(r)
<Residue 'ALA', 99>
>>> r.name = "UNK"
>>> print(r)
<Residue 'UNK', 99>
>>> rnew = u.select_atoms("resid 99").residues[0]
>>> print(rnew)
<Residue 'UNK', 99>

will typically work as expected. When working with collections such as AtomGroup or ResidueGroup it is generally better to use properties such as AtomGroup.resnames or ResidueGroup.resnames to modify their items’ attributes.

There are two cases when it is very important to use these collective properties:

Because residues are determined by the Atom.resid and segments by Atom.segid, the above properties take extra care to rebuild the list of segments and residues. Alternatively, the same effect can be obtained using the corresponding setter method, e.g. AtomGroup.set_resids().

Note

Setting any of AtomGroup.resids, ResidueGroup.resids, AtomGroup.segids, ResidueGroup.segids can change the topology: they can split or merge residues or segments.

Splitting/merging of residues is probably not very useful because no chemical rearrangements are carried out. Manipulating segments might be more useful in order to add additional structure to a Universe and provide instant segment selectors for interactive work:

u.select_atoms("protein").set_segids("protein")
u.select_atoms("resname POPE or resname POPC").set_segids("lipids")
u.select_atoms("resname SOL").set_segids("water")
u.select_atoms("resname NA or resname CL").set_segids("ions")

u.protein.n_residues
water_oxygens = u.water.OW

The setter methods have the additional advantage that they can assign lists. For instance, many MD codes number residues consecutively starting from 1. However, the original structure might be missing a few atoms at the N-terminus. Let’s say that the first residue is really residue 10. In order to store the canonical residue IDs (“resnum”) one could the use

protein = u.select_atoms("protein").residues
protein.set_resnums(protein.resnums + 9)

One can then use protein.select("resnum 42") to select the residue that has the canonical residue id 42 (instead of resid 33).

One can also read the resids directly from an original PDB file:

orig = MDAnalysis.Universe("2jln.pdb")
protein.set_resnums(orig.select_atoms("protein").resids)

7.2.3. Working with Topologies

If the topology file given to the Universe had bonding information then this will have been loaded into the Universe as Universe.bonds Universe.angles Universe.dihedrals and Universe.impropers.

If your topology file does not have this information, it is still possible to construct it based on the positions of the atoms and assumed vdw radii for these atoms. See MDAnalysis.AtomGroup.guess_bonds() and MDAnalysis.topology.core.guess_bonds() for details.

This Topology information is stored in MDAnalysis.core.topologyobjects.TopologyGroup objects. These are designed to be analogous to the AtomGroup container for Atoms.

For examples working with a box of ethanol:

>>> import MDAnalysis as mda
>>> u = mda.Universe('ethanol.gro', guess_bonds=True)
>>> u.bonds
<TopologyGroup containing 11784 Bonds>
>>> u.bonds.types()  # view available types
[('O', 'H'), ('C', 'O'), ('C', 'H'), ('C', 'C')]
>>> u.bonds.select_bonds(('C', 'O'))  # return all C-O bonds from the group
<TopologyGroup containing 1473 Bonds>

Bonds are categorised based on the types of atoms. This is done in a way so that type (a, b, c) is equivalent to (c, b, a) ie. bonds are reversible. For example:

>>> u.angles.types()
[('C', 'C', 'H'),
 ('H', 'C', 'H'),
 ('C', 'O', 'H'),
 ('C', 'C', 'O'),
 ('H', 'C', 'O')]

There are only C-C-H bonds and no H-C-C bonds. Selection however is aware that sometimes types are reversed:

>>> u.angles.select_bonds(('H', 'C', 'C'))  # note reversal of type
<TopologyGroup containing 7365 Angles>

TopologyGroups can be combined and indexed:

>>> tg = u.angles.select_bonds(('C', 'C', 'O')) + u.angles.select_bonds(('C', 'O', 'H'))
>>> tg.types()
[('C', 'O', 'H'), ('C', 'C', 'O')]
>>> tg[:100]
<TopologyGroup containing 100 Angles>

Finally, TopologyGroups are linked to some fast Cython calculation methods to determine bond lengths and angle sizes:

>>> tg.values()
array([ 1.88042373,  1.95928987,  1.74770012, ...,  1.79306789,
        1.95522678,  1.88881045])

7.2.4. Combining objects: system building

It is often convenient to combined multiple groups of atoms into a single object. If they are contained in a single Universe then the methods described above (especially manipulating the segments) might be useful. However, if the atoms reside in different universes, the Merge() function can be used.

7.2.4.1. Merging

In the following example for Merge(), protein, ligand, and solvent were externally prepared in three different PDB files. They are loaded into separate Universe objects (where they could be further manipulated, e.g. renumbered, relabeled, rotated, ...) The Merge() command is used to combine all of them together:

import MDAnalysis
u1 = MDAnalysis.Universe("protein.pdb")
u2 = MDAnalysis.Universe("ligand.pdb")
u3 = MDAnalysis.Universe("solvent.pdb")
u = MDAnalysis.Merge(u1.select_atoms("protein"), u2.atoms, u3.atoms)
u.atoms.write("system.pdb")

The complete system is then written out to a new PDB file.

7.2.4.2. Replicating

It is also possible to replicate a molecule to build a system with multiple copies of the same molecule. In the example, we replicate an AdK molecule and then translate and rotate the second copy:

import MDAnalysis; from MDAnalysis.tests.datafiles import *
u = MDAnalysis.Universe(PSF, DCD)
p = u.select_atoms("protein")
m = MDAnalysis.Merge(p,p)

# now renumber resids and segids for each copy

# first copy of the protein (need to use atom indices because currently that's the only reliable property in the
merged universe)
p1 = m.select_atoms("bynum 1:3341")
# second copy
p2 = m.select_atoms("bynum 3342:6682")

p1.set_segid("A")
p2.set_segid("B")
p2.residues.set_resid(p2.residues.resids + p1.residues.resids[-1])  # increment resids for p2 with the last
resid from p1

# you must regenerate the selections after modifying them (see notes in the docs!)
# because the changed resids are not reflected in the selection (due to how residues are referenced internally)
p1 = m.select_atoms("segid A")       # or as before:  m.select_atoms("bynum 1:3341")
p2 = m.select_atoms("segid B")

# rotate and translate
p2.rotateby(180, [0,0,1])
p2.translate([50,0,0])

Note that we have to manually set the residue numbers (resids) and segment identifies because Merge() simply concatenates the existing atoms and only ensures that all data structures are contained in the new merged universe.

7.2.5. Classes and functions

class MDAnalysis.core.AtomGroup.Universe(*args, **kwargs)[source]

The MDAnalysis Universe contains all the information describing the system.

The system always requires a topology file — in the simplest case just a list of atoms. This can be a CHARMM/NAMD PSF file or a simple coordinate file with atom informations such as XYZ, PDB, Gromacs GRO, or CHARMM CRD. See Table of Supported Topology Formats for what kind of topologies can be read.

A trajectory provides coordinates; the coordinates have to be ordered in the same way as the list of atoms in the topology. A trajectory can be a single frame such as a PDB, CRD, or GRO file, or it can be a MD trajectory (in CHARMM/NAMD/LAMMPS DCD, Gromacs XTC/TRR, or generic XYZ format). See Table of supported coordinate formats for what can be read as a “trajectory”.

As a special case, when the topology is a file that contains atom information and coordinates (such as XYZ, PDB, GRO or CRD, see Table of supported coordinate formats) then the coordinates are immediately loaded from the “topology” file unless a trajectory is supplied.

Examples for setting up a universe:

u = Universe(topology, trajectory)          # read system from file(s)
u = Universe(pdbfile)                       # read atoms and coordinates from PDB or GRO
u = Universe(topology, [traj1, traj2, ...]) # read from a list of trajectories
u = Universe(topology, traj1, traj2, ...)   # read from multiple trajectories

Load new data into a universe (replaces old trajectory and does not append):

u.load_new(trajectory)                      # read from a new trajectory file

Select atoms, with syntax similar to CHARMM (see select_atoms for details):

u.select_atoms(...)

Attributes:

Note

If atom attributes such as element, mass, or charge are not explicitly provided in the topology file then MDAnalysis tries to guess them (see MDAnalysis.topology.tables). This does not always work and if you require correct values (e.g. because you want to calculate the center of mass) then you need to make sure that MDAnalysis gets all the information needed.

Changed in version 0.7.5: Can also read multi-frame PDB files with the PrimitivePDBReader.

Changed in version 0.8: Parse arbitrary number of arguments as a single topology file and a sequence of trajectories.

Changed in version 0.9.0: Topology information now loaded lazily, but can be forced with build_topology(). Changed bonds attribute to be a TopologyGroup. Added angles and torsions attribute as TopologyGroup. Added fragments to Universe cache

Changed in version 0.11.0: make_anchor(), remove_anchor(), is_anchor, and anchor_name were added to support the pickling/unpickling of AtomGroup. Deprecated selectAtoms() in favour of select_atoms().

Initialize the central MDAnalysis Universe object.

Arguments:
topologyfile

A CHARMM/XPLOR PSF topology file, PDB file or Gromacs GRO file; used to define the list of atoms. If the file includes bond information, partial charges, atom masses, ... then these data will be available to MDAnalysis. A “structure” file (PSF, PDB or GRO, in the sense of a topology) is always required.

coordinatefile

A trajectory (such as CHARMM DCD, Gromacs XTC/TRR/GRO, XYZ, XYZ.bz2) or a PDB that will provide coordinates, possibly multiple frames. If a list of filenames is provided then they are sequentially read and appear as one single trajectory to the Universe. The list can contain different file formats.

Deprecated since version 0.8: Do not use the coordinatefile keyword argument, just provide trajectories as positional arguments.

permissive

currently only relevant for PDB files: Set to True in order to ignore most errors and read typical MD simulation PDB files; set to False to read with the Bio.PDB reader, which can be useful for real Protein Databank PDB files. None selects the MDAnalysis default (which is set in MDAnalysis.core.flags) [None]

topology_format

provide the file format of the topology file; None guesses it from the file extension [None] Can also pass a subclass of MDAnalysis.topology.base.TopologyReader to define a custom reader to be used on the topology file.

format

provide the file format of the coordinate or trajectory file; None guesses it from the file extension. Note that this keyword has no effect if a list of file names is supplied because the “chained” reader has to guess the file format for each individual list member. [None] Can also pass a subclass of MDAnalysis.coordinates.base.Reader to define a custom reader to be used on the trajectory file.

guess_bonds

Once Universe has been loaded, attempt to guess the connectivity between atoms. This will populate the .bonds .angles and .dihedrals attributes of the Universe.

vdwradii

For use with guess_bonds. Supply a dict giving a vdwradii for each atom type which are used in guessing bonds.

is_anchor

When unpickling instances of MDAnalysis.core.AtomGroup.AtomGroup existing Universes are searched for one where to anchor those atoms. Set to False to prevent this Universe from being considered. [True]

anchor_name

Setting to other than None will cause MDAnalysis.core.AtomGroup.AtomGroup instances pickled from the Universe to only unpickle if a compatible Universe with matching anchor_name is found. is_anchor will be ignored in this case but will still be honored when unpickling MDAnalysis.core.AtomGroup.AtomGroup instances pickled with *anchor_name*==``None``. [None]

This class tries to do the right thing:

  1. If file with topology and coordinate information (such as PDB, GRO, CRD, ...) is provided instead of a topology file and no coordinatefile then the coordinates are taken from the first file. Thus you can load a functional universe with

    u = Universe('1ake.pdb')
    

    If you want to specify the coordinate file format yourself you can do so using the format keyword:

    u = Universe('1ake.ent1', format='pdb')
    
  2. If only a topology file without coordinate information is provided one will have to load coordinates manually using Universe.load_new(). The file format of the topology file can be explicitly set with the topology_format keyword.

Changed in version 0.7.4: New topology_format and format parameters to override the file format detection.

Changed in version 0.10.0: Added 'guess_bonds' keyword to cause topology to be guessed on Universe creation. Deprecated 'bonds' keyword, use 'guess_bonds' instead.

Changed in version 0.11.0: Added the is_anchor and anchor_name keywords for finer behavior control when unpickling instances of MDAnalysis.core.AtomGroup.AtomGroup.

angles

Returns a TopologyGroup of all angles in the Universe

New in version 0.9.0.

Changed in version 0.9.2: Now can return empty TopologyGroup

bonds

Returns a TopologyGroup of all bonds in the Universe.

Changed in version 0.9.0: Now a lazily built TopologyGroup

Changed in version 0.9.2: Now can return empty TopologyGroup

build_topology()[source]

Bond angle and dihedral information is lazily constructed into the Universe.

This method forces all this information to be loaded.

coord

Reference to current timestep and coordinates of universe.

The raw trajectory coordinates are Universe.coord.positions, represented as a numpy.float32 array.

Because coord is a reference to a Timestep, it changes its contents while one is stepping through the trajectory.

Note

In order to access the coordinates it is probably better to use the AtomGroup.coordinates() method; for instance, all coordinates of the Universe as a numpy array: Universe.atoms.coordinates().

dihedrals

Returns a TopologyGroup of all dihedrals in the Universe

New in version 0.9.0.

Changed in version 0.9.2: Now can return empty TopologyGroup

dimensions

Current dimensions of the unitcell

fragments

Read only tuple of fragments in the Universe

impropers

Returns a TopologyGroup of all improper dihedrals in the Universe

New in version 0.9.0.

Changed in version 0.9.2: Now can return empty TopologyGroup

is_anchor

Whether this Universe will be checked for anchoring when unpickling MDAnalysis.core.AtomGroup.AtomGroup instances

load_new(filename, **kwargs)[source]

Load coordinates from filename, using the suffix to detect file format.

Arguments:
filename

the coordinate file (single frame or trajectory) or a list of filenames, which are read one after another.

permissive

currently only relevant for PDB files: Set to True in order to ignore most errors and read typical MD simulation PDB files; set to False to read with the Bio.PDB reader, which can be useful for real Protein Databank PDB files. None selects the MDAnalysis default (which is set in MDAnalysis.core.flags) [None]

format

provide the file format of the coordinate or trajectory file; None guesses it from the file extension. Note that this keyword has no effect if a list of file names is supplied because the “chained” reader has to guess the file format for each individual list member [None] Can also pass a subclass of MDAnalysis.coordinates.base.Reader to define a custom reader to be used on the trajectory file.

kwargs

Other kwargs are passed to the trajectory reader (only for advanced use)

Returns:

(filename, trajectory_format) or None if filename == None

Raises:

TypeError if trajectory format can not be determined or no appropriate trajectory reader found

Changed in version 0.8: If a list or sequence that is provided for filename only contains a single entry then it is treated as single coordinate file. This has the consequence that it is not read by the ChainReader but directly by its specialized file format reader, which typically has more features than the ChainReader.

make_anchor()[source]

Add this Universe to the list where anchors are searched for when unpickling MDAnalysis.core.AtomGroup.AtomGroup instances. Silently proceeds if it is already on the list.

remove_anchor()[source]

Remove this Universe from the list where anchors are searched for when unpickling MDAnalysis.core.AtomGroup.AtomGroup instances. Silently proceeds if it is already not on the list.

selectAtoms(*args, **kwds)

selectAtoms is deprecated, use select_atoms instead!

Select atoms using a CHARMM selection string.

Returns an AtomGroup with atoms sorted according to their index in the psf (this is to ensure that there aren’t any duplicates, which can happen with complicated selections).

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

Subselections can be grouped with parentheses.

Example::
>>> sel = universe.select_atoms("segid DMPC and not ( name H* or name O* )")
>>> sel
<AtomGroup with 3420 atoms>
>>> universe.select_atoms("around 10 group notHO", notHO=sel)
<AtomGroup with 1250 atoms>

Note

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.

See also

Selection Commands for further details and examples.

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.
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.

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 or resname 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.

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

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 MDAnalysis.Universe are consecutively numbered, and the index runs from 1 up to the total number of atoms.

Preexisting selections

group group-name
selects the atoms in the AtomGroup passed to the function as an argument named group-name. Only the atoms common to group-name and the instance select_atoms() was called from will be considered. 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 they were defined in an appropriate Universe.
fullgroup group-name
just like the group keyword with the difference that all the atoms of group-name are included. The resulting selection may therefore have atoms that were initially absent from the instance select_atoms() was called from.

Changed in version 0.7.4: Added resnum selection.

Changed in version 0.8.1: Added group and fullgroup selections.

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

Select atoms using a CHARMM selection string.

Returns an AtomGroup with atoms sorted according to their index in the psf (this is to ensure that there aren’t any duplicates, which can happen with complicated selections).

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

Subselections can be grouped with parentheses.

Example::
>>> sel = universe.select_atoms("segid DMPC and not ( name H* or name O* )")
>>> sel
<AtomGroup with 3420 atoms>
>>> universe.select_atoms("around 10 group notHO", notHO=sel)
<AtomGroup with 1250 atoms>

Note

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.

See also

Selection Commands for further details and examples.

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.
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.

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 or resname 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.

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

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 MDAnalysis.Universe are consecutively numbered, and the index runs from 1 up to the total number of atoms.

Preexisting selections

group group-name
selects the atoms in the AtomGroup passed to the function as an argument named group-name. Only the atoms common to group-name and the instance select_atoms() was called from will be considered. 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 they were defined in an appropriate Universe.
fullgroup group-name
just like the group keyword with the difference that all the atoms of group-name are included. The resulting selection may therefore have atoms that were initially absent from the instance select_atoms() was called from.

Changed in version 0.7.4: Added resnum selection.

Changed in version 0.8.1: Added group and fullgroup selections.

trajectory

Reference to trajectory reader object containing trajectory data.

MDAnalysis.core.AtomGroup.Merge(*args)[source]

Return a Universe from two or more AtomGroup instances.

AtomGroup instances can come from different Universes, or come directly from a select_atoms() call.

It can also be used with a single AtomGroup if the user wants to, for example, re-order the atoms in the Universe.

Arguments:One or more AtomGroup instances.
Returns:an instance of Universe
Raises:ValueError for too few arguments or if an AtomGroup is empty and TypeError if arguments are not AtomGroup instances.

Example

In this example, protein, ligand, and solvent were externally prepared in three different PDB files. They are loaded into separate Universe objects (where they could be further manipulated, e.g. renumbered, relabeled, rotated, ...) The Merge() command is used to combine all of them together:

u1 = Universe("protein.pdb")
u2 = Universe("ligand.pdb")
u3 = Universe("solvent.pdb")
u = Merge(u1.select_atoms("protein"), u2.atoms, u3.atoms)
u.atoms.write("system.pdb")

The complete system is then written out to a new PDB file.

Note

Merging does not create a full trajectory but only a single structure even if the input consists of one or more trajectories.

class MDAnalysis.core.AtomGroup.AtomGroup(atoms)[source]

A group of atoms.

ag = universe.select_atoms(atom-list)

The AtomGroup contains a list of atoms; typically, a AtomGroup is generated from a selection. It is build from any list-like collection of Atom instances. It is also possible to create an empty AtomGroup from an empty list.

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 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 containing (ag[0], ag[3], ag[-1], ag[10], ag[3]).

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).

Atoms can also be accessed in a Pythonic fashion by using the atom name as an attribute. For instance,

ag.CA

will provide a AtomGroup of all CA atoms in the group. These instant selector attributes are auto-generated for each atom name encountered in the group.

Note

The name-attribute instant selector access to atoms is mainly meant for quick interactive work. Thus it either returns a single Atom if there is only one matching atom, or a new AtomGroup for multiple matches. This makes it difficult to use the feature consistently in scripts but it is much better for interactive work.

AtomGroup instances are bound to a Universe, but only through the weak reference Atom has. The connection is lost as soon as the Universe goes out of scope.

AtomGroups can also be pickled and unpickled. When pickling, the Atom indices are serialized alongside the Universe number-of-atoms, and topology and trajectory filenames. When unpickling, all built Universes will be searched for one with matching number-of-atoms and filenames. Finer control over which Universe gets chosen to base the unpickling on can be exerted using the is_anchor and anchor_name flags upon Universe creation or the methods/attributes Universe.make_anchor(), Universe.remove_anchor(), and Universe.anchor_name.

References for analysis methods

[Dima2004](1, 2, 3) Dima, R. I., & Thirumalai, D. (2004). Asymmetry in the shapes of folded and denatured states of proteins. J Phys Chem B, 108(21), 6564-6570. doi:10.1021/jp037128y

See also

Universe

Changed in version 0.7.6: An empty AtomGroup can be created and no longer raises a NoDataError.

Changed in version 0.9.0: The size at which cache is used for atom lookup is now stored as variable _atomcache_size within the class. Added fragments manged property. Is a lazily built, cached entry, similar to residues.

Changed in version 0.11.0: AtomGroups can now be pickled and unpickled provided compatible Universes are available. The follow methods were changed to properties: indices, masses, charges, names, types, radii, resids, resnames, resnums, segids. Added altLocs and serials properties and setters. Torsions and torsion renamed to dihedral. The bond, angle, dihedral and improper methods were removed and replaced with properties of the same name which return the corresponding object. Deprecated selectAtoms in favour of select_atoms. Setters are now plural to match property names. Properties referring to residue (resids, resnames, resnums) or segment [segids] properties now yield arrays of length equal to n_atoms

_atoms

immutable list of references to the atoms Atom in the group

_rebuild_caches()[source]

Rebuild all AtomGroup caches.

A number of lists and attributes are cached. These caches are lazily built the first time they are needed. When editing the topology it might happen that not all caches were synced properly (even though that this is supposed to happen eventually). In this case the user can manually force a complete cache rebuild.

Currently the following caches are used:

See also

_clear_caches()

New in version 0.7.5.

Changed in version 0.9.0: Added bonds/angles/torsions/impropers to rebuild. Reworked how things are rebuilt to avoid code duplication.

_clear_caches(*args)[source]

Clear cache for all args.

If no args are provided, all caches are cleared.

New in version 0.8.

_fill_cache(name, value)[source]

Populate _cache[name] with value.

New in version 0.8.

align_principalAxis(*args, **kwds)

align_principalAxis is deprecated, use align_principal_axis instead!

Align principal axis with index axis with vector.

Arguments:
axis

Index of the principal axis (0, 1, or 2), as produced by principal_axes().

vector

A 3D vector such as the z-axis ([0,0,1]); can be anything that looks like a list with three entries.

To align the long axis of a channel (the first principal axis, i.e. axis = 0) with the z-axis:

u.atoms.align_principal_axis(0, [0,0,1])
u.atoms.write("aligned.pdb")
align_principal_axis(axis, vector)[source]

Align principal axis with index axis with vector.

Arguments:
axis

Index of the principal axis (0, 1, or 2), as produced by principal_axes().

vector

A 3D vector such as the z-axis ([0,0,1]); can be anything that looks like a list with three entries.

To align the long axis of a channel (the first principal axis, i.e. axis = 0) with the z-axis:

u.atoms.align_principal_axis(0, [0,0,1])
u.atoms.write("aligned.pdb")
altLocs

numpy array of the altLocs for all atoms in this group

New in version 0.11.0.

angle

This AtomGroup represented as an Angle object

Returns:A MDAnalysis.core.topologyobjects.Angle object
Raises:ValueError if the AtomGroup is not length 3

New in version 0.11.0.

angles

All the angles in this AtomGroup

Note that these angles might extend out of the AtomGroup, to select only angles which are entirely contained by the AtomGroup use u.angles.atomgroup_intersection(ag, strict=True)

New in version 0.9.0.

Changed in version 0.10.0: Now always returns a (possibly empty) TopologyGroup

asphericity(**kwargs)[source]

Asphericity.

See [Dima2004] for background information.

Keywords:
pbc

True: Move all atoms within primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

New in version 0.7.7.

Changed in version 0.8: Added pbc keyword

atoms

AtomGroup of all atoms in this group.

If this is a AtomGroup then it returns itself. Otherwise, it will return a new AtomGroup based on all Atom instances contained.

Apply :func:`list to atoms or use _atoms if you really only need a list of individual Atom instances.

bbox(**kwargs)[source]

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]
Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Returns:[[xmin, ymin, zmin], [xmax, ymax, zmax]]

New in version 0.7.2.

Changed in version 0.8: Added pbc keyword

bfactors

Crystallographic B-factors (from PDB) in A**2.

bond

This AtomGroup represented as a Bond object

Returns:A MDAnalysis.core.topologyobjects.Bond object
Raises:ValueError if the AtomGroup is not length 2

New in version 0.11.0.

bonds

All the bonds in this AtomGroup

Note that these bonds might extend out of the AtomGroup, to select only bonds which are entirely contained by the AtomGroup use u.bonds.atomgroup_intersection(ag, strict=True)

New in version 0.9.0.

Changed in version 0.10.0: Now always returns a (possibly empty) TopologyGroup

bsphere(**kwargs)[source]

Return the bounding sphere of the selection.

The sphere is calculated relative to the centre of geometry.

Keywords:
pbc

True: Move all atoms within primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Returns:(R, [xcen,ycen,zcen])

New in version 0.7.3.

Changed in version 0.8: Added pbc keyword

centerOfGeometry(*args, **kwds)

centerOfGeometry is deprecated, use center_of_geometry instead!

Center of geometry (also known as centroid) of the selection.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

centerOfMass(*args, **kwds)

centerOfMass is deprecated, use center_of_mass instead!

Center of mass of the selection.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

center_of_geometry(**kwargs)[source]

Center of geometry (also known as centroid) of the selection.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

center_of_mass(**kwargs)[source]

Center of mass of the selection.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

centroid(**kwargs)

Center of geometry (also known as centroid) of the selection.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

charges

Array of partial charges of the atoms (as defined in the topology)

Changed in version 0.11.0: Now a property

coordinates(ts=None, copy=False, dtype=<type 'numpy.float32'>)

NumPy array of the coordinates.

See also

positions and get_positions()

Deprecated since version 0.7.6: In new scripts use AtomGroup.get_positions() preferrably.

dihedral

This AtomGroup represented as a Dihedral object

Returns:A MDAnalysis.core.topologyobjects.Dihedral object
Raises:ValueError if the AtomGroup is not length 4

New in version 0.11.0.

dihedrals

All the dihedrals in this AtomGroup

Note that these dihedrals might extend out of the AtomGroup, to select only dihedrals which are entirely contained by the AtomGroup use u.dihedrals.atomgroup_intersection(ag, strict=True)

New in version 0.9.0.

Changed in version 0.10.0: Now always returns a (possibly empty) TopologyGroup

dimensions

Dimensions of the Universe to which the group belongs, at the current time step.

forces

Forces on the atoms in the AtomGroup.

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

For more control use the get_forces() and set_forces() methods.

New in version 0.7.7.

fragments

Read-only list of fragments.

Contains all fragments that any Atom in this AtomGroup is part of, the contents of the fragments may extend beyond the contents of this AtomGroup.

get_forces(ts=None, copy=False, dtype=<type 'numpy.float32'>)[source]

Get a NumPy array of the atomic forces (if available). Currently only supported for Gromacs .trr trajectories.

Keywords:
ts

If ts is provided then positions are read from that Timestep instead of the one from the current trajectory belonging to this universe. The ts is indexed with the indices returned by indices() and it is the user’s responsibility to provide a time step that has the appropriate dimensions.

copy

True: always make a copy (slow), False: Try to return a array view or reference (faster); note that for passing coordinates to C-code it can be necessary to use a copy [False]

dtype

NumPy Data type of the array; the default is usually entirely appropriate. Most C-code actually requires the default [numpy.float32]

Forces can also be directly obtained from the attribute forces.

Forces can be directly set with set_forces() or by assigning to forces.

New in version 0.7.7.

get_positions(ts=None, copy=False, dtype=<type 'numpy.float32'>)[source]

Get a NumPy array of the coordinates.

Keywords:
ts

If ts is provided then positions are read from that Timestep instead of the one from the current trajectory belonging to this universe. The ts is indexed with the indices returned by indices() and it is the user’s responsibility to provide a time step that has the appropriate dimensions.

copy

True: always make a copy (slow), False: Try to return a array view or reference (faster); note that for passing coordinates to C-code it can be necessary to use a copy [False]

dtype

NumPy Data type of the array; the default is usually entirely appropriate. Most C-code actually requires the default [numpy.float32]

Coordinates can also be directly obtained from the attribute positions.

Coordinates can be directly set with set_positions() or by assigning to positions.

This method is identical with coordinates() but named differently for symmetry with with set_positions().

New in version 0.7.6.

get_velocities(ts=None, copy=False, dtype=<type 'numpy.float32'>)[source]

NumPy array of the velocities.

Raises a NoDataError if the underlying Timestep does not contain velocities.

See also AtomGroup.set_velocities() and attribute access through AtomGroup.velocities.

New in version 0.7.6.

guess_bonds(vdwradii=None)[source]

Guess all the bonds that exist within this AtomGroup and add to Universe.

Keywords:
vdwradii

Pass a dict relating atom types to vdwradii.

New in version 0.10.0.

improper

This AtomGroup represented as an ImproperDihedral object

Returns:A MDAnalysis.core.topologyobjects.ImproperDihedral object
Raises:ValueError if the AtomGroup is not length 4

New in version 0.11.0.

impropers

All the improper dihedrals in this AtomGroup

Note that these improper dihedrals might extend out of the AtomGroup, to select only dihedrals which are entirely contained by the AtomGroup use u.impropers.atomgroup_intersection(ag, strict=True)

New in version 0.9.0.

Changed in version 0.10.0: Now always returns a (possibly empty) TopologyGroup

indices

Array of all Atom.index in the group.

These indices are 0-based and can be used to directly index Universe.atoms or the coordinate array MDAnalysis.coordinates.base.Timestep.positions.

Note

This property is read only

Changed in version 0.11.0: Now a property

masses

Array of atomic masses (as defined in the topology)

Changed in version 0.11.0: Now a property

momentOfInertia(*args, **kwds)

momentOfInertia is deprecated, use moment_of_inertia instead!

Tensor of inertia as 3x3 NumPy array.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

moment_of_inertia(**kwargs)[source]

Tensor of inertia as 3x3 NumPy array.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

n_atoms

Total number of atoms in the group

n_residues

Total number of residues in the group

n_segments

Total number of segments in the group

names

Returns an array of atom names.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property

packIntoBox(*args, **kwds)

packIntoBox is deprecated, use pack_into_box instead!

Shift all atoms in this group to be within the primary unit cell.

AtomGroup.packintobox([box, [inplace=True]])

Keywords:
box

Unit cell to move atoms inside of.

inplace

True: Change coordinates in place and return False: Only return the coordinates

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.

By default the coordinates are changed in place and returned

New in version 0.8.

pack_into_box(box=None, inplace=True)[source]

Shift all atoms in this group to be within the primary unit cell.

AtomGroup.packintobox([box, [inplace=True]])

Keywords:
box

Unit cell to move atoms inside of.

inplace

True: Change coordinates in place and return False: Only return the coordinates

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.

By default the coordinates are changed in place and returned

New in version 0.8.

positions

Coordinates of the atoms in the AtomGroup.

The positions can be changed by assigning an array of the appropriate shape, i.e. either Nx3 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 particles to the origin).

For more control use the get_positions() and set_positions() methods.

New in version 0.7.6.

principalAxes(*args, **kwds)

principalAxes is deprecated, use principal_axes instead!

Calculate the principal axes from the moment of inertia.

e1,e2,e3 = AtomGroup.principal_axes()

The eigenvectors are sorted by eigenvalue, i.e. the first one corresponds to the highest eigenvalue and is thus the first principal axes.

Keywords:
pbc

True: Move all atoms within primary unit cell before calculation

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Returns:numpy.array v with v[0] as first, v[1] as second, and v[2] as third eigenvector.

Changed in version 0.8: Added pbc keyword

principal_axes(**kwargs)[source]

Calculate the principal axes from the moment of inertia.

e1,e2,e3 = AtomGroup.principal_axes()

The eigenvectors are sorted by eigenvalue, i.e. the first one corresponds to the highest eigenvalue and is thus the first principal axes.

Keywords:
pbc

True: Move all atoms within primary unit cell before calculation

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Returns:numpy.array v with v[0] as first, v[1] as second, and v[2] as third eigenvector.

Changed in version 0.8: Added pbc keyword

radii

Array of atomic radii (as defined in the PQR file)

Changed in version 0.11.0: Now a property

radiusOfGyration(*args, **kwds)

radiusOfGyration is deprecated, use radius_of_gyration instead!

Radius of gyration.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

radius_of_gyration(**kwargs)[source]

Radius of gyration.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

Changed in version 0.8: Added pbc keyword

resids

Returns an array of residue numbers.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

residues

Read-only list of Residue objects.

A ResidueGroup of all residues that contain atoms in this group.

Changed in version 0.9.0: Now returns strictly a ResidueGroup of the unique Residues that Atoms in this group belong to.

resnames

Returns an array of residue names.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

resnums

Returns an array of canonical residue numbers.

New in version 0.7.4.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

rotate(R)[source]

Apply a rotation matrix R to the selection’s coordinates.

AtomGroup.rotate(R)

\(\mathsf{R}\) is a 3x3 orthogonal matrix that transforms a vector \(\mathbf{x} \rightarrow \mathbf{x}'\):

\[\mathbf{x}' = \mathsf{R}\mathbf{x}\]
rotateby(angle, axis, point=None)[source]

Apply a rotation to the selection’s coordinates.

AtomGroup.rotateby(angle,axis[,point])

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}\).

Arguments:
angle

rotation angle in degrees

axis

rotation axis vector, a 3-tuple, list, or array, or a 2-tuple of two MDAnalysis objects from which the axis is calculated as the vector from the first to the second center of geometry.

point

point on the rotation axis; by default (None) the center of geometry of the selection is chosen, or, if axis is a tuple of selections, it defaults to the first point of the axis. point can be a 3-tuple, list, or array or a MDAnalysis object (in which case its centroid() is used).

Returns:

The 4x4 matrix which consists of the rotation matrix M[:3,:3] and the translation vector M[:3,3].

segids

Returns an array of segment names.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

segments

Read-only list of Segment objects.

A SegmentGroup of all segments that contain atoms in this group.

Changed in version 0.9.0: Now strictly returns a SegmentGroup of a set of the Segment instances from this AtomGroup

selectAtoms(*args, **kwds)

selectAtoms is deprecated, use select_atoms instead!

Selection of atoms using the MDAnalysis selection syntax.

AtomGroup.select_atoms(selection[,selection[,...]], [groupname=atomgroup[,groupname=atomgroup[,...]]])

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

Selection of atoms using the MDAnalysis selection syntax.

AtomGroup.select_atoms(selection[,selection[,...]], [groupname=atomgroup[,groupname=atomgroup[,...]]])

sequence(**kwargs)[source]

Returns the amino acid sequence.

The format of the sequence is selected with the keyword format:

format description
‘SeqRecord’ Bio.SeqRecord.SeqRecord (default)
‘Seq’ Bio.Seq.Seq
‘string’ string

The sequence is returned by default (keyword format = 'SeqRecord') as a Bio.SeqRecord.SeqRecord instance, which can then be further processed. In this case, all keyword arguments (such as the id string or the name or the description) are directly passed to Bio.SeqRecord.SeqRecord.

If the keyword format is set to 'Seq', all kwargs are ignored and a Bio.Seq.Seq instance is returned. The difference to the record is that the record also contains metadata and can be directly used as an input for other functions in Bio.

If the keyword format is set to 'string', all kwargs are ignored and a Python string is returned.

Example: Write FASTA file

Use Bio.SeqIO.write(), which takes sequence records:

import Bio.SeqIO

# get the sequence record of a protein component of a Universe
protein = u.select_atoms("protein")
record = protein.sequence(id="myseq1", name="myprotein")

Bio.SeqIO.write(record, "single.fasta", "fasta")

A FASTA file with multiple entries can be written with

Bio.SeqIO.write([record1, record2, ...], "multi.fasta", "fasta")
Keywords:

format

  • "string": return sequence as a string of 1-letter codes
  • "Seq": return a Bio.Seq.Seq instance
  • "SeqRecord": return a Bio.SeqRecord.SeqRecord instance

Default is "SeqRecord"

id

Sequence ID for SeqRecord (should be different for different sequences)

name

Name of the protein.

description

Short description of the sequence.

kwargs

Any other keyword arguments that are understood by Bio.SeqRecord.SeqRecord.

Raises:

ValueError if a residue name cannot be converted to a 1-letter IUPAC protein amino acid code; make sure to only select protein residues. Raises TypeError if an unknown format is selected.

New in version 0.9.0.

serials

numpy array of the serials for all atoms in this group

New in version 0.11.0.

set(name, value, **kwargs)

Set attribute name to value for all atoms in the AtomGroup.

If value is a sequence of the same length as the AtomGroup then each Atom‘s property name is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms to distinct values by providing a sequence or iterable.

set_altLoc(*args, **kwds)

set_altLoc is deprecated, use set_altLocs instead!

Set the altLocs to altLoc for **all atoms* in the AtomGroup.

If altLoc is a sequence of the same length as the AtomGroup then each Atom.altLoc is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.11.0.

set_altLocs(altLoc)[source]

Set the altLocs to altLoc for **all atoms* in the AtomGroup.

If altLoc is a sequence of the same length as the AtomGroup then each Atom.altLoc is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.11.0.

set_bfactor(*args, **kwds)

set_bfactor is deprecated, use set_bfactors instead!

Set the atom bfactors to float bfactor for all atoms in the AtomGroup.

If bfactor is a sequence of the same length as the AtomGroup then each Atom.bfactor is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_bfactors(bfactor)[source]

Set the atom bfactors to float bfactor for all atoms in the AtomGroup.

If bfactor is a sequence of the same length as the AtomGroup then each Atom.bfactor is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_charge(*args, **kwds)

set_charge is deprecated, use set_charges instead!

Set the partial charges to float charge for all atoms in the AtomGroup.

If charge is a sequence of the same length as the AtomGroup then each Atom.charge is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_charges(charge)[source]

Set the partial charges to float charge for all atoms in the AtomGroup.

If charge is a sequence of the same length as the AtomGroup then each Atom.charge is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_forces(forces, ts=None)[source]

Set the forces for all atoms in the group.

Arguments:
forces

a Nx3 NumPy numpy.ndarray where N is the number of atoms in this atom group.

Keywords:
ts

Timestep, defaults to None and then the current time step is used.

Note

If the group contains N atoms and force is a single vector (i.e. an array of length 3) then all N atom positions are set to force (due to NumPy’s broadcasting rules), as described for forces.

See also get_forces() and attribute access through forces.

New in version 0.7.7.

set_mass(*args, **kwds)

set_mass is deprecated, use set_masses instead!

Set the atom masses to float mass for all atoms in the AtomGroup.

If mass is a sequence of the same length as the AtomGroup then each Atom.mass is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_masses(mass)[source]

Set the atom masses to float mass for all atoms in the AtomGroup.

If mass is a sequence of the same length as the AtomGroup then each Atom.mass is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_name(*args, **kwds)

set_name is deprecated, use set_names instead!

Set the atom names to string for all atoms in the AtomGroup.

If value is a sequence of the same length as the AtomGroup then each Atom.name is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_names(name)[source]

Set the atom names to string for all atoms in the AtomGroup.

If value is a sequence of the same length as the AtomGroup then each Atom.name is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_positions(coords, ts=None)[source]

Set the positions for all atoms in the group.

Arguments:
coords

a Nx3 NumPy numpy.ndarray where N is the number of atoms in this atom group.

Keywords:
ts

Timestep, defaults to None and then the current time step is used.

Note

If the group contains N atoms and coord is a single point (i.e. an array of length 3) then all N atom positions are set to coord (due to NumPy’s broadcasting rules), as described for positions.

See also get_positions() and attribute access through positions.

New in version 0.7.6.

set_radii(radius)[source]

Set the atom radii to float radius for all atoms in the AtomGroup.

If radius is a sequence of the same length as the AtomGroup then each Atom.radius is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_radius(*args, **kwds)

set_radius is deprecated, use set_radii instead!

Set the atom radii to float radius for all atoms in the AtomGroup.

If radius is a sequence of the same length as the AtomGroup then each Atom.radius is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resid(*args, **kwds)

set_resid is deprecated, use set_resids instead!

Set the resids to integer resid for all atoms in the AtomGroup.

If resid is a sequence of the same length as the AtomGroup then each Atom.resid is set to the corresponding value together with the Residue.id of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resids can change the topology.

Assigning the same resid to multiple residues will merge these residues. Assigning different resid to atoms in the same residue will split a residue (and potentially merge with another one).

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable and can change the topology via MDAnalysis.topology.core.build_residues().

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resids(resid)[source]

Set the resids to integer resid for all atoms in the AtomGroup.

If resid is a sequence of the same length as the AtomGroup then each Atom.resid is set to the corresponding value together with the Residue.id of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resids can change the topology.

Assigning the same resid to multiple residues will merge these residues. Assigning different resid to atoms in the same residue will split a residue (and potentially merge with another one).

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable and can change the topology via MDAnalysis.topology.core.build_residues().

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resname(*args, **kwds)

set_resname is deprecated, use set_resnames instead!

Set the resnames to string resname for all atoms in the AtomGroup.

If resname is a sequence of the same length as the AtomGroup then each Atom.resname is set to the corresponding value together with the Residue.name of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resnames(resname)[source]

Set the resnames to string resname for all atoms in the AtomGroup.

If resname is a sequence of the same length as the AtomGroup then each Atom.resname is set to the corresponding value together with the Residue.name of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resnum(*args, **kwds)

set_resnum is deprecated, use set_resnums instead!

Set the resnums to resnum for all atoms in the AtomGroup.

If resnum is a sequence of the same length as the AtomGroup then each Atom.resnum is set to the corresponding value together with the Residue.resnum of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resnum will not affect the topology: you can have multiple residues with the same resnum.

See also

set_resid()

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resnums(resnum)[source]

Set the resnums to resnum for all atoms in the AtomGroup.

If resnum is a sequence of the same length as the AtomGroup then each Atom.resnum is set to the corresponding value together with the Residue.resnum of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resnum will not affect the topology: you can have multiple residues with the same resnum.

See also

set_resid()

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_segid(*args, **kwds)

set_segid is deprecated, use set_segids instead!

Set the segids to segid for all atoms in the AtomGroup.

If segid is a sequence of the same length as the AtomGroup then each Atom.segid is set to the corresponding value together with the Segment.id of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

set_segid() can change the topology.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Stale caches are problematic; though it can be expensive, changing segid results in Segments being regenerated

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_segids(segid)[source]

Set the segids to segid for all atoms in the AtomGroup.

If segid is a sequence of the same length as the AtomGroup then each Atom.segid is set to the corresponding value together with the Segment.id of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

set_segid() can change the topology.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Stale caches are problematic; though it can be expensive, changing segid results in Segments being regenerated

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_serial(*args, **kwds)

set_serial is deprecated, use set_serials instead!

Set the serials to serial for all atoms in the AtomGroup.

If serial is a sequence of the same length as the AtomGroup then each Atom.serial is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.11.0.

set_serials(serial)[source]

Set the serials to serial for all atoms in the AtomGroup.

If serial is a sequence of the same length as the AtomGroup then each Atom.serial is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.11.0.

set_type(*args, **kwds)

set_type is deprecated, use set_types instead!

Set the atom types to atype for all atoms in the AtomGroup.

If atype is a sequence of the same length as the AtomGroup then each Atom.atype is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_types(atype)[source]

Set the atom types to atype for all atoms in the AtomGroup.

If atype is a sequence of the same length as the AtomGroup then each Atom.atype is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_velocities(v, ts=None)[source]

Assign the velocities v to the timestep.

Raises a NoDataError if the underlying Timestep does not contain velocities.

See also AtomGroup.get_velocities() and AtomGroup.velocities for attribute access.

New in version 0.7.6.

shapeParameter(*args, **kwds)

shapeParameter is deprecated, use shape_parameter instead!

Shape parameter.

See [Dima2004] for background information.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

New in version 0.7.7.

Changed in version 0.8: Added pbc keyword

shape_parameter(**kwargs)[source]

Shape parameter.

See [Dima2004] for background information.

Keywords:
pbc

True: Move all atoms within the primary unit cell before calculation [False]

Note

The MDAnalysis.core.flags flag use_pbc when set to True allows the pbc flag to be used by default.

New in version 0.7.7.

Changed in version 0.8: Added pbc keyword

split(level)[source]

Split atomgroup into a list of atomgroups by level.

level can be “atom”, “residue”, “segment”. .. versionadded:: 0.9.0

totalCharge(*args, **kwds)

totalCharge is deprecated, use total_charge instead!

Sum of all partial charges (must be defined in topology).

totalMass(*args, **kwds)

totalMass is deprecated, use total_mass instead!

Total mass of the selection (masses are taken from the topology or guessed).

total_charge()[source]

Sum of all partial charges (must be defined in topology).

total_mass()[source]

Total mass of the selection (masses are taken from the topology or guessed).

transform(M)[source]

Apply homogenous transformation matrix M to the coordinates.

The matrix M must be a 4x4 matrix, with the rotation in R = `M[:3,:3]` and the translation in t = M[:3,3].

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

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

Apply translation vector t to the selection’s coordinates.

>>> AtomGroup.translate(t)
>>> AtomGroup.translate((A, B))

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

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

The translation can also be given as a tuple of two MDAnalysis objects such as two selections (selA, selB), i.e. two AtomGroup, or two Atom. The translation vector is computed as the difference between the centers of geometry (centroid) of B and A:

t = B.centroid() - A.centroid()
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.

types

Returns an array of atom types.

New in version 0.9.0.

Changed in version 0.11.0: Now a property

universe

The universe to which the atoms belong (read-only).

velocities

NumPy array of the velocities of the atoms in the group.

If the trajectory does not contain velocity information then a NoDataError is raised.

New in version 0.7.5.

Deprecated since version 0.7.6: In 0.8 this will become an attribute! You can already use get_velocities() and set_velocities().

Changed in version 0.8: Became an attribute.

wrap(compound='atoms', center='com', box=None)[source]

Shift the contents of this AtomGroup back into the unit cell.

This is a more powerful version of pack_into_box(), allowing groups of atoms to be kept together through the process.

Keywords:
compound

The group which will be kept together through the shifting process. [atoms]

Possible options:

  • atoms
  • group - This AtomGroup
  • residues
  • segments
  • fragments
center

How to define the center of a given group of atoms [com]

box

Unit cell information. If not provided, the values from Timestep will be used.

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 all atoms from the compound are not inside the unit cell, but rather the center of the compound is. Compounds available for use are atoms, residues, segments and fragments

center allows the definition of the center of each group to be specified. 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, an the dimensions information from the current Timestep will be used.

Note

wrap with all default keywords is identical to pack_into_box()

New in version 0.9.2.

write(filename=None, format='PDB', filenamefmt='%(trjname)s_%(frame)d', **kwargs)[source]

Write AtomGroup to a file.

AtomGroup.write(filename[,format])

Keywords:
filename

None: create TRJNAME_FRAME.FORMAT from filenamefmt [None]

format

PDB, CRD, GRO, VMD (tcl), PyMol (pml), Gromacs (ndx) CHARMM (str) Jmol (spt); case-insensitive and can also be supplied as the filename extension [PDB]

filenamefmt

format string for default filename; use substitution tokens ‘trjname’ and ‘frame’ [“%(trjname)s_%(frame)d”]

bonds

how to handle bond information, especially relevant for PDBs; default is "conect".

  • "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

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

write_selection(filename=None, format='vmd', filenamefmt='%(trjname)s_%(frame)d', **kwargs)[source]

Write AtomGroup selection to a file to be used in another programme.

Keywords:
filename

None: create TRJNAME_FRAME.FORMAT from filenamefmt

format

output file format: VMD (tcl), PyMol (pml), Gromacs (ndx), CHARMM (str); can also be supplied as the filename extension. Case insensitive. [vmd]

filenamefmt

format string for default filename; use ‘%(trjname)s’ and ‘%(frame)s’ placeholders; the extension is set according to the format [“%(trjname)s_%(frame)d”]

kwargs

additional keywords are passed on to the appropriate SelectionWriter

Deprecated since version 0.9.0: Use write()

class MDAnalysis.core.AtomGroup.Atom(index, name, type, resname, resid, segid, mass, charge, residue=None, segment=None, radius=None, bfactor=None, resnum=None, serial=None, altLoc=None, universe=None)[source]

A class representing a single atom.

Atom is the basic building block of all larger data structures in MDAnalysis, in particular of the AtomGroup.

An Atom is typically generated by a topology reader from MDAnalysis.topology.

For performance reasons, only a predefined number of attributes are included (and thus it is not possible to add attributes “on the fly”; they have to be included in the class definition).

An Atom is bound to a particular Universe, but via a weak reference only. This means that the Atom, and any AtomGroup it is in, are only relevant while the originating Universe is in scope.

index

atom index

segid

name of the segment

resid

residue number

resnum

canonical residue number as, for instance, used in the original PDB file

New in version 0.7.4.

resname

residue name

residue

Residue object containing the atoms

id

atom number inside the residue

name

string, short name

type

string or number (from force field), describing the atom type; stored as a string.

Changed in version 0.7.6: The Atom.type attribute is always stored as a string.

mass

float, in atomic mass units (u)

charge

float, in electron charges (e)

radius

Born-radius for electrostatic calculations. (Only if read from a PQR file with the PQRReader.)

altLoc

Alternate location indicator (as used in ATOM_ records in PDB files)

New in version 0.9.0.

bonds

A MDAnalysis.core.topologyobjects.TopologyGroup of all Bond instances that contains all the bonds that this atom participates in.

New in version 0.8.1.

angles

A MDAnalysis.core.topologyobjects.TopologyGroup of all Angle instances that contains all the angles that this atom participates in.

New in version 0.8.1.

dihedrals

A MDAnalysis.core.topologyobjects.TopologyGroup of all Dihedral instances that contains all the dihedrals that this atom participates in.

New in version 0.8.1.

impropers

A MDAnalysis.core.topologyobjects.TopologyGroup of all MDAnalysis.core.topologyobjects.ImproperDihedral instances that this atom is present in.

New in version 0.8.1.

angles

A TopologyGroup of the angles that this Atom is in

Changed in version 0.9.0: Changed to managed property

bonded_atoms

An AtomGroup of the Atoms that this Atom is bonded to.

New in version 0.9.

bonds

A TopologyGroup of the bonds that this Atom is in

Changed in version 0.9.0: Changed to managed property

centroid()[source]

The centroid of an atom is its position, Atom.position.

dihedrals

A TopologyGroup of the dihedrals that this Atom is in

Changed in version 0.9.0: Changed to managed property

Changed in version 0.11.0: Renamed to dihedrals (was torsions)

force

Current force of the atom.

Returns:a (3,) shape numpy array

A NoDataError is raised if the trajectory does not contain velocities.

New in version 0.9.2.

fragment

The fragment that this Atom is part of

New in version 0.9.0.

impropers

A TopologyGroup of the improper dihedrals that this Atom is in

Changed in version 0.9.0: Changed to managed property

number

The index of this atom

pos

coordinates of the atom

Get the current cartesian coordinates of the atom (read-only).

Deprecated since version 0.8: use position

position

coordinates of the atom

Get the current cartesian coordinates of the atom.

Returns:a (3,) shape numpy array
universe

a pointer back to the Universe

velocity

Current velocity of the atom.

Returns:a (3,) shape numpy array

A NoDataError is raised if the trajectory does not contain velocities.

New in version 0.7.5.

class MDAnalysis.core.AtomGroup.Residue(name, id, atoms, resnum=None)[source]

A group of atoms corresponding to a residue.

Pythonic access to atoms:
  • Using a atom name as attribute returns the matching atom (a Atom instance), i.e. r.name. Example:

    >>> from MDAnalysis.tests.datafiles import PSF,DCD
    >>> u = Universe(PSF,DCD)
    >>> print(u.s4AKE.r1.CA)  # C-alpha of M1
    < Atom 5: name 'CA' of type '22' of resname 'MET', resid 1 and segid '4AKE'>
    
  • r['name'] or r[id] - returns the atom corresponding to that name

Data:
Residue.name

Three letter residue name.

Residue.id

Numeric (integer) resid, taken from the topology.

Residue.resnum

Numeric canonical residue id (e.g. as used in the PDB structure).

Note

Creating a Residue modifies the underlying Atom instances. Each Atom can only belong to a single Residue.

Changed in version 0.7.4: Added Residue.resnum attribute and resnum keyword argument.

chi1_selection()[source]

AtomGroup corresponding to the chi1 sidechain dihedral N-CA-CB-CG.

Returns:4-atom selection in the correct order. If no CB and/or CG is found then this method returns None.

New in version 0.7.5.

omega_selection()[source]

AtomGroup corresponding to the omega protein backbone dihedral CA-C-N’-CA’.

omega describes the -C-N- peptide bond. Typically, it is trans (180 degrees) although cis-bonds (0 degrees) are also occasionally observed (especially near Proline).

Returns:4-atom selection in the correct order. If no C’ found in the previous residue (by resid) then this method returns None.
phi_selection()[source]

AtomGroup corresponding to the phi protein backbone dihedral C’-N-CA-C.

Returns:4-atom selection in the correct order. If no C’ found in the previous residue (by resid) then this method returns None.
psi_selection()[source]

AtomGroup corresponding to the psi protein backbone dihedral N-CA-C-N’.

Returns:4-atom selection in the correct order. If no N’ found in the following residue (by resid) then this method returns None.
class MDAnalysis.core.AtomGroup.ResidueGroup(residues)[source]

A group of residues.

Pythonic access to atoms:
  • Using a atom name as attribute returns a list of all atoms (a AtomGroup) of the same name. Example:

    >>> from MDAnalysis.tests.datafiles import PSF,DCD
    >>> u = Universe(PSF,DCD)
    >>> print(u.s4AKE.MET.CA)  # C-alpha of all Met
    <AtomGroup with 6 atoms>
    
Data:ResidueGroup._residues

Initialize the ResidueGroup with a list of Residue instances.

resids

Returns an array of residue numbers.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

resnames

Returns an array of residue names.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

resnums

Returns an array of canonical residue numbers.

New in version 0.7.4.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

segids

Returns an array of segment names.

Note

uses the segment of the first atom in each residue for the segment name returned

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

set(name, value, **kwargs)

Set attribute name to value for all residues in the ResidueGroup.

If value is a sequence of the same length as the ResidueGroup (AtomGroup.residues) then each Residue‘s property name is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the ResidueGroup then a ValueError is raised.

New in version 0.7.5.

Changed in version 0.8: Can set residues to distinct values by providing a sequence or iterable.

set_resid(*args, **kwds)

set_resid is deprecated, use set_resids instead!

Set the resids to integer resid for all residues in the

ResidueGroup.

If resid is a sequence of the same length as the ResidueGroup then each Atom.resid is set to the corresponding value together with the Residue.id of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resids can change the topology.

Assigning the same resid to multiple residues will merge these residues. The new residue name will be the name of the first old residue in the merged residue.

Warning

The values of this ResidueGroup are not being changed. You must create a new ResidueGroup from the Universe — only Atom instances are changed, everything else is derived from these atoms.

New in version 0.8.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resids(resid)[source]

Set the resids to integer resid for all residues in the ResidueGroup.

If resid is a sequence of the same length as the ResidueGroup then each Atom.resid is set to the corresponding value together with the Residue.id of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resids can change the topology.

Assigning the same resid to multiple residues will merge these residues. The new residue name will be the name of the first old residue in the merged residue.

Warning

The values of this ResidueGroup are not being changed. You must create a new ResidueGroup from the Universe — only Atom instances are changed, everything else is derived from these atoms.

New in version 0.8.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resname(*args, **kwds)

set_resname is deprecated, use set_resnames instead!

Set the resnames to string resname for all residues in the

ResidueGroup.

If resname is a sequence of the same length as the ResidueGroup then each Atom.resname is set to the corresponding value together with the Residue.name of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resnames(resname)[source]

Set the resnames to string resname for all residues in the ResidueGroup.

If resname is a sequence of the same length as the ResidueGroup then each Atom.resname is set to the corresponding value together with the Residue.name of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resnum(*args, **kwds)

set_resnum is deprecated, use set_resnums instead!

Set the resnums to resnum for all residues in the ResidueGroup.

If resnum is a sequence of the same length as the ResidueGroup then each Atom.resnum is set to the corresponding value together with the Residue.resnum of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resnum will not affect the topology: you can have multiple residues with the same resnum.

See also

set_resid()

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_resnums(resnum)[source]

Set the resnums to resnum for all residues in the ResidueGroup.

If resnum is a sequence of the same length as the ResidueGroup then each Atom.resnum is set to the corresponding value together with the Residue.resnum of the residue the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

Changing resnum will not affect the topology: you can have multiple residues with the same resnum.

See also

set_resid()

New in version 0.7.4.

Changed in version 0.7.5: Also changes the residues.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

class MDAnalysis.core.AtomGroup.Segment(name, residues)[source]

A group of residues corresponding to one segment of the topology.

Pythonic access to residues:

  • The attribute rN returns the N-th residue Residue of the segment (numbering starts at N=1). Example:

    >>> from MDAnalysis.tests.datafiles import PSF,DCD
    >>> u = Universe(PSF,DCD)
    >>> print(u.s4AKE.r1)
    <Residue 'MET', 1>
    
  • Using a residue name as attribute returns a list of all residues (a ResidueGroup) of the same name. Example:

    >>> from MDAnalysis.tests.datafiles import PSF,DCD
    >>> u = Universe(PSF,DCD)
    >>> print(u.s4AKE.CYS)
    <ResidueGroup [<Residue 'CYS', 77>]>
    >>> print(u.s4AKE.MET)
    <ResidueGroup [<Residue 'MET', 1>, <Residue 'MET', 21>, <Residue 'MET', 34>, <Residue 'MET', 53>,
    <Residue 'MET', 96>, <Residue 'MET', 174>]>
    
Data:Segment.name is the segid from the topology or the chain identifier when loaded from a PDB

Initialize a Segment with segid name from a list of Residue instances.

id

Segment id (alias for Segment.name)

class MDAnalysis.core.AtomGroup.SegmentGroup(segments)[source]

A group of segments.

Pythonic access to segments:
  • Using a segid as attribute returns the segment. Because of python language rule, any segid starting with a non-letter character is prefixed with ‘s’, thus ‘4AKE’ –> ‘s4AKE’.

    Example:

    >>> from MDAnalysis.tests.datafiles import PSF,DCD
    >>> u = Universe(PSF,DCD)
    >>> print(u.atoms.segments.s4AKE)  # segment 4AKE
    <AtomGroup with 3314 atoms>
    
  • Indexing the group returns the appropriate segment.

Data:SegmentGroup._segments

Initialize the SegmentGroup with a list of Segment instances.

segids

Returns an array of segment names.

Changed in version 0.8: Returns a numpy.ndarray

Changed in version 0.11.0: Now a property and returns array of length len(self)

set(name, value, **kwargs)

Set attribute name to value for all Segment in this AtomGroup.

If value is a sequence of the same length as the SegmentGroup (AtomGroup.residues) then each Segment‘s property name is set to the corresponding value. If value is neither of length 1 (or a scalar) nor of the length of the SegmentGroup then a ValueError is raised.

New in version 0.8.

set_segid(*args, **kwds)

set_segid is deprecated, use set_segids instead!

Set the segids to segid for all atoms in the SegmentGroup.

If segid is a sequence of the same length as the SegmentGroup then each Atom.segid is set to the corresponding value together with the Segment.id of the segment the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

set_segid() can change the topology.

This can be used to join segments or to break groups into multiple disjoint segments. Note that each Atom can only belong to a single Segment.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

set_segids(segid)[source]

Set the segids to segid for all atoms in the SegmentGroup.

If segid is a sequence of the same length as the SegmentGroup then each Atom.segid is set to the corresponding value together with the Segment.id of the segment the atom belongs to. If value is neither of length 1 (or a scalar) nor of the length of the AtomGroup then a ValueError is raised.

Note

set_segid() can change the topology.

This can be used to join segments or to break groups into multiple disjoint segments. Note that each Atom can only belong to a single Segment.

New in version 0.7.4.

Changed in version 0.8: Can set atoms and residues to distinct values by providing a sequence or iterable.

Changed in version 0.11.0: Made plural to make consistent with corresponding property

MDAnalysis.core.AtomGroup.as_Universe(*args, **kwargs)[source]

Return a universe from the input arguments.

  1. If the first argument is a universe, just return it:

    as_Universe(universe) --> universe
    
  2. Otherwise try to build a universe from the first or the first and second argument:

    as_Universe(PDB, **kwargs) --> Universe(PDB, **kwargs)
    as_Universe(PSF, DCD, **kwargs) --> Universe(PSF, DCD, **kwargs)
    as_Universe(*args, **kwargs) --> Universe(*args, **kwargs)
    
Returns:an instance of Universe
exception MDAnalysis.core.AtomGroup.SelectionError[source]

Raised when a atom selection failed.

exception MDAnalysis.core.AtomGroup.SelectionWarning[source]

Warning indicating a possible problem with a selection.

exception MDAnalysis.core.AtomGroup.NoDataError[source]

Raised when empty input is not allowed or required data are missing.