11.2.1. Core Topology object — MDAnalysis.core.topology

New in version 0.16.0.

Topology is the core object that holds all topology information.

TODO: Add in-depth discussion.

Notes

For developers: In MDAnalysis 0.16.0 this new topology system was introduced and discussed as issue #363; this issue contains key information and discussions on the new system. The issue number 363 is also being used as a short-hand in discussions to refer to the new topology system.

11.2.1.1. Classes

class MDAnalysis.core.topology.Topology(n_atoms=1, n_res=1, n_seg=1, attrs=None, atom_resindex=None, residue_segindex=None)[source]

In-memory, array-based topology database.

The topology model of MDanalysis features atoms, which must each be a member of one residue. Each residue, in turn, must be a member of one segment. The details of maintaining this heirarchy, and mappings of atoms to residues, residues to segments, and vice-versa, are handled internally by this object.

Parameters:
  • n_atoms (int) – number of atoms in topology. Must be larger then 1 at each level
  • n_residues (int) – number of residues in topology. Must be larger then 1 at each level
  • n_segments (int) – number of segments in topology. Must be larger then 1 at each level
  • attrs (TopologyAttr objects) – components of the topology to be included
  • atom_resindex (array) – 1-D array giving the resindex of each atom in the system
  • residue_segindex (array) – 1-D array giving the segindex of each residue in the system
add_Residue(segment, **new_attrs)[source]
Returns:
Return type:residx of the new Residue
Raises:NoDataError – If not all data was provided. This error is raised before any
add_TopologyAttr(topologyattr)[source]

Add a new TopologyAttr to the Topology.

Parameters:topologyattr (TopologyAttr) –
copy()[source]

Return a deepcopy of this Topology

guessed_attributes

A list of the guessed attributes in this topology

read_attributes

A list of the attributes read from the topology

class MDAnalysis.core.topology.TransTable(n_atoms, n_residues, n_segments, atom_resindex=None, residue_segindex=None)[source]

Membership tables with methods to translate indices across levels.

There are three levels; Atom, Residue and Segment. Each Atom must belong in a Residue, each Residue must belong to a Segment.

When translating upwards, eg finding which Segment a Residue belongs in, a single numpy array is returned. When translating downwards, two options are available; a concatenated result (suffix _1) or a list for each parent object (suffix _2d).

Parameters:
  • n_atoms (int) – number of atoms in topology
  • n_residues (int) – number of residues in topology
  • n_segments (int) – number of segments in topology
  • atom_resindex (1-D array) – resindex for each atom in the topology; the number of unique values in this array must be <= n_residues, and the array must be length n_atoms; giving None defaults to placing all atoms in residue 0
  • residue_segindex (1-D array) – segindex for each residue in the topology; the number of unique values in this array must be <= n_segments, and the array must be length n_residues; giving None defaults to placing all residues in segment 0
n_atoms

number of atoms in topology

Type:int
n_residues

number of residues in topology

Type:int
n_segments

number of segments in topology

Type:int
size

tuple describing the shape of the TransTable

atoms2residues(aix)[source]

Returns the residue index for many atom indices

residues2atoms_1d(rix)[source]

All atoms in the residues represented by rix

residues2atoms_2d(rix)[source]

List of atom indices for each residue in rix

residues2segments(rix)[source]

Segment indices for each residue in rix

segments2residues_1d(six)[source]

Similar to residues2atoms_1d

segments2residues_2d(six)[source]

Similar to residues2atoms_2d

atoms2segments(aix)[source]

Segment indices for each atom in aix

segments2atoms_1d(six)[source]

Similar to residues2atoms_1d

segments2atoms_2d(six)[source]

Similar to residues2atoms_2d

atoms2residues(aix)[source]

Get residue indices for each atom.

Parameters:aix (array) – atom indices
Returns:rix – residue index for each atom
Return type:array
atoms2segments(aix)[source]

Get segment indices for each atom.

Parameters:aix (array) – atom indices
Returns:rix – segment index for each atom
Return type:array
copy()[source]

Return a deepcopy of this Transtable

move_atom(aix, rix)[source]

Move aix to be in rix

move_residue(rix, six)[source]

Move rix to be in six

residues2atoms_1d(rix)[source]

Get atom indices collectively represented by given residue indices.

Parameters:rix (array) – residue indices
Returns:aix – indices of atoms present in residues, collectively
Return type:array
residues2atoms_2d(rix)[source]

Get atom indices represented by each residue index.

Parameters:rix (array) – residue indices
Returns:raix – each element corresponds to a residue index, in order given in rix, with each element being an array of the atom indices present in that residue
Return type:list
residues2segments(rix)[source]

Get segment indices for each residue.

Parameters:rix (array) – residue indices
Returns:six – segment index for each residue
Return type:array
segments2atoms_1d(six)[source]

Get atom indices collectively represented by given segment indices.

Parameters:six (array) – segment indices
Returns:aix – sorted indices of atoms present in segments, collectively
Return type:array
segments2atoms_2d(six)[source]

Get atom indices represented by each segment index.

Parameters:six (array) – residue indices
Returns:saix – each element corresponds to a segment index, in order given in six, with each element being an array of the atom indices present in that segment
Return type:list
segments2residues_1d(six)[source]

Get residue indices collectively represented by given segment indices

Parameters:six (array) – segment indices
Returns:rix – sorted indices of residues present in segments, collectively
Return type:array
segments2residues_2d(six)[source]

Get residue indices represented by each segment index.

Parameters:six (array) – residue indices
Returns:srix – each element corresponds to a segment index, in order given in six, with each element being an array of the residue indices present in that segment
Return type:list
size

The shape of the table, (n_atoms, n_residues, n_segments)

11.2.1.2. Helper functions

MDAnalysis.core.topology.make_downshift_arrays(upshift, nparents)[source]

From an upwards translation table, create the opposite direction

Turns a many to one mapping (eg atoms to residues) to a one to many mapping (residues to atoms)

Parameters:
  • upshift (array_like) – Array of integers describing which parent each item belongs to
  • nparents (integer) – Total number of parents that exist.
Returns:

downshift – An array of arrays, each containing the indices of the children of each parent. Length nparents + 1

Return type:

array_like (dtype object)

Examples

To find the residue to atom mappings for a given atom to residue mapping:

>>> atom2res = np.array([0, 1, 0, 2, 2, 0, 2])
>>> make_downshift_arrays(atom2res)
array([array([0, 2, 5]), array([1]), array([3, 4, 6]), None], dtype=object)

Entry 0 corresponds to residue 0 and says that this contains atoms 0, 2 & 5

Notes

The final entry in the return array will be None to ensure that the dtype of the array is object.

Warning

This means negative indexing should never be used with these arrays.