mdtraj.Topology¶
-
class
mdtraj.
Topology
¶ Topology stores the topological information about a system.
The structure of a Topology object is similar to that of a PDB file. It consists of a set of Chains (often but not always corresponding to polymer chains). Each Chain contains a set of Residues, and each Residue contains a set of Atoms. In addition, the Topology stores a list of which atom pairs are bonded to each other.
Atom and residue names should follow the PDB 3.0 nomenclature for all molecules for which one exists.
Examples
>>> topology = md.load('example.pdb').topology >>> print(topology) <mdtraj.Topology with 1 chains, 3 residues, 22 atoms, 21 bonds at 0x105a98e90> >>> table, bonds = topology.to_dataframe() >>> print(table.head()) serial name element resSeq resName chainID 0 0 H1 H 1 CYS 0 1 1 CH3 C 1 CYS 0 2 2 H2 H 1 CYS 0 3 3 H3 H 1 CYS 0 4 4 C C 1 CYS 0 >>> # rename residue "CYS" to "CYSS" >>> table[table['residue'] == 'CYS']['residue'] = 'CYSS' >>> print(table.head()) serial name element resSeq resName chainID 0 0 H1 H 1 CYSS 0 1 1 CH3 C 1 CYSS 0 2 2 H2 H 1 CYSS 0 3 3 H3 H 1 CYSS 0 4 4 C C 1 CYSS 0 >>> t2 = md.Topology.from_dataframe(table, bonds)
- Attributes
Methods
add_atom
(self, name, element, residue[, serial])Create a new Atom and add it to the Topology.
add_bond
(self, atom1, atom2[, type, order])Create a new bond and add it to the Topology.
add_chain
(self)Create a new Chain and add it to the Topology.
add_residue
(self, name, chain[, resSeq, …])Create a new Residue and add it to the Topology.
atom
(self, index)Get a specific atom by index.
atoms_by_name
(self, name)Iterator over all Atoms in the Topology with a specified name
chain
(self, index)Get a specific chain by index.
copy
(self)Return a copy of the topology
create_disulfide_bonds
(self, positions)Identify disulfide bonds based on proximity and add them to the Topology.
create_standard_bonds
(self)Create bonds based on the atom and residue names for all standard residue types.
delete_atom_by_index
(self, index)Delete an Atom from the topology.
find_molecules
(self)Identify molecules based on bonds.
from_dataframe
(atoms[, bonds])Create a mdtraj topology from a pandas data frame
from_openmm
(value)Create a mdtraj topology from an OpenMM topology
guess_anchor_molecules
(self)Guess anchor molecules for imaging
insert_atom
(self, name, element, residue[, …])Create a new Atom and insert it into the Topology at a specific position.
join
(self, other)Join two topologies together
residue
(self, index)Get a specific residue by index.
select
(self, selection_string)Execute a selection against the topology
select_atom_indices
(self[, selection])Get the indices of biologically-relevant groups by name.
select_expression
(self, selection_string)Translate a atom selection expression into a pure python expression.
select_pairs
(self[, selection1, selection2])Generate unique pairs of atom indices.
subset
(self, atom_indices)Create a new Topology from a subset of the atoms in an existing topology.
to_bondgraph
(self)Create a NetworkX graph from the atoms and bonds in this topology
to_dataframe
(self)Convert this topology into a pandas dataframe
to_fasta
(self[, chain])Convert this topology into FASTA string
to_openmm
(self[, traj])Convert this topology into OpenMM topology
-
__init__
(self)¶ Create a new Topology object
Methods
__init__
(self)Create a new Topology object
add_atom
(self, name, element, residue[, serial])Create a new Atom and add it to the Topology.
add_bond
(self, atom1, atom2[, type, order])Create a new bond and add it to the Topology.
add_chain
(self)Create a new Chain and add it to the Topology.
add_residue
(self, name, chain[, resSeq, …])Create a new Residue and add it to the Topology.
atom
(self, index)Get a specific atom by index.
atoms_by_name
(self, name)Iterator over all Atoms in the Topology with a specified name
chain
(self, index)Get a specific chain by index.
copy
(self)Return a copy of the topology
create_disulfide_bonds
(self, positions)Identify disulfide bonds based on proximity and add them to the Topology.
create_standard_bonds
(self)Create bonds based on the atom and residue names for all standard residue types.
delete_atom_by_index
(self, index)Delete an Atom from the topology.
find_molecules
(self)Identify molecules based on bonds.
from_dataframe
(atoms[, bonds])Create a mdtraj topology from a pandas data frame
from_openmm
(value)Create a mdtraj topology from an OpenMM topology
guess_anchor_molecules
(self)Guess anchor molecules for imaging
insert_atom
(self, name, element, residue[, …])Create a new Atom and insert it into the Topology at a specific position.
join
(self, other)Join two topologies together
residue
(self, index)Get a specific residue by index.
select
(self, selection_string)Execute a selection against the topology
select_atom_indices
(self[, selection])Get the indices of biologically-relevant groups by name.
select_expression
(self, selection_string)Translate a atom selection expression into a pure python expression.
select_pairs
(self[, selection1, selection2])Generate unique pairs of atom indices.
subset
(self, atom_indices)Create a new Topology from a subset of the atoms in an existing topology.
to_bondgraph
(self)Create a NetworkX graph from the atoms and bonds in this topology
to_dataframe
(self)Convert this topology into a pandas dataframe
to_fasta
(self[, chain])Convert this topology into FASTA string
to_openmm
(self[, traj])Convert this topology into OpenMM topology
Attributes
Iterator over all Atoms in the Topology.
Iterator over all bonds (each represented as a tuple of two Atoms) in the Topology.
Iterator over all Chains in the Topology.
Get the number of atoms in the Topology
Get the number of bonds in the Topology
Get the number of chains in the Topology
Get the number of residues in the Topology.
Iterator over all Residues in the Topology.
-
add_atom
(self, name, element, residue, serial=None)¶ Create a new Atom and add it to the Topology.
- Parameters
- namestr
The name of the atom to add
- elementmdtraj.element.Element
The element of the atom to add
- residuemdtraj.topology.Residue
The Residue to add it to
- serialint
Serial number associated with the atom.
- Returns
- atommdtraj.topology.Atom
the newly created Atom
-
add_bond
(self, atom1, atom2, type=None, order=None)¶ Create a new bond and add it to the Topology.
- Parameters
- atom1mdtraj.topology.Atom
The first Atom connected by the bond
- atom2mdtraj.topology.Atom
The second Atom connected by the bond
- typemdtraj.topology.Singleton or None, Default: None, Optional
Bond type of the bond, or None if not known/provided
- order1, 2, 3 or None, Default: None, Optional
Characteristic order of the bond if known, defaults None
-
add_chain
(self)¶ Create a new Chain and add it to the Topology.
- Returns
- chainmdtraj.topology.Chain
the newly created Chain
-
add_residue
(self, name, chain, resSeq=None, segment_id='')¶ Create a new Residue and add it to the Topology.
- Parameters
- namestr
The name of the residue to add
- chainmdtraj.topology.Chain
The Chain to add it to
- resSeqint, optional
Residue sequence number, such as from a PDB record. These sequence numbers are arbitrary, and do not necessarily start at 0 (or 1). If not supplied, the resSeq attribute will be set to the residue’s sequential (0 based) index.
- segment_idstr, optional
A label for the segment to which this residue belongs
- Returns
- residuemdtraj.topology.Residue
The newly created Residue
-
atom
(self, index)¶ Get a specific atom by index. These indices start from zero.
- Parameters
- indexint
The index of the atom to select.
- Returns
- atomAtom
The index-th atom in the topology.
-
property
atoms
¶ Iterator over all Atoms in the Topology.
- Returns
- atomitergenerator
Iterator over all Atoms in the Topology.
-
atoms_by_name
(self, name)¶ Iterator over all Atoms in the Topology with a specified name
- Parameters
- namestr
The particular atom name of interest.
- Returns
- atomitergenerator
Examples
>>> for atom in topology.atoms_by_name('CA'): ... print(atom)
-
property
bonds
¶ Iterator over all bonds (each represented as a tuple of two Atoms) in the Topology.
- Returns
- bonditergenerator
Iterator over all bonds between Atoms in the Trajectory. Each bond can then be iterated to get the atoms. I.e.: for atom1, atom2 in Topology.bonds yields iterator over pairs of Atoms
for bond in Topology.bonds yields iterator over bonds
-
chain
(self, index)¶ Get a specific chain by index. These indices start from zero.
- Parameters
- indexint
The index of the chain to select.
- Returns
- chainChain
The index-th chain in the topology.
-
property
chains
¶ Iterator over all Chains in the Topology.
- Returns
- chainiterlistiterator
Iterator over all Chains in the Topology.
-
copy
(self)¶ Return a copy of the topology
- Returns
- outTopology
A copy of this topology
-
create_disulfide_bonds
(self, positions)¶ Identify disulfide bonds based on proximity and add them to the Topology.
- Parameters
- positionslist
The list of atomic positions based on which to identify bonded atoms
-
create_standard_bonds
(self)¶ Create bonds based on the atom and residue names for all standard residue types.
-
delete_atom_by_index
(self, index)¶ Delete an Atom from the topology.
- Parameters
- indexint
The index of the atom to be removed.
-
find_molecules
(self)¶ Identify molecules based on bonds.
A molecule is defined as a set of atoms that are connected to each other by bonds. This method uses the list of bonds to divide up the Topology’s atoms into molecules.
- Returns
- moleculeslist of sets
Each entry represents one molecule, and is the set of all Atoms in that molecule
-
classmethod
from_dataframe
(atoms, bonds=None)¶ Create a mdtraj topology from a pandas data frame
- Parameters
- atomspandas.DataFrame
The atoms in the topology, represented as a data frame. This data frame should have columns “serial” (atom index), “name” (atom name), “element” (atom’s element), “resSeq” (index of the residue) “resName” (name of the residue), “chainID” (index of the chain), and optionally “segmentID”, following the same conventions as wwPDB 3.0 format.
- bondsnp.ndarray, shape=(n_bonds, 4) or (n_bonds, 2), dtype=float, Optional
The bonds in the topology, represented as a n_bonds x 4 or n_bonds x 2 size array of the indices of the atoms involved, type, and order of each bond, represented as floats. Type and order are optional. Specifying bonds here is optional. To create standard protein bonds, you can use create_standard_bonds to “fill in” the bonds on your newly created Topology object, although type and order of bond are not computed if that method is used.
See also
-
classmethod
from_openmm
(value)¶ Create a mdtraj topology from an OpenMM topology
- Parameters
- valuesimtk.openmm.app.Topology
An OpenMM topology that you wish to convert to a mdtraj topology.
-
guess_anchor_molecules
(self)¶ Guess anchor molecules for imaging
- Returns
- anchor_moleculeslist of atom sets
List of anchor molecules
See also
-
insert_atom
(self, name, element, residue, index=None, rindex=None, serial=None)¶ Create a new Atom and insert it into the Topology at a specific position.
- Parameters
- namestr
The name of the atom to add
- elementmdtraj.element.Element
The element of the atom to add
- residuemdtraj.topology.Residue
The Residue to add it to
- indexint
If provided, the desired index for this atom within the topology. Existing atoms with indices >= index will be pushed back.
- rindexint
If provided, the desired position for this atom within the residue
- serialint
Serial number associated with the atom. This has nothing to do with the actual ordering and is solely for PDB labeling purposes.
- Returns
- atommdtraj.topology.Atom
the newly created Atom
-
join
(self, other)¶ Join two topologies together
- Parameters
- otherTopology
Another topology object
- Returns
- outTopology
A joint topology, with all of the atoms/residues/chains/bonds in each of the individual topologies
-
property
n_atoms
¶ Get the number of atoms in the Topology
-
property
n_bonds
¶ Get the number of bonds in the Topology
-
property
n_chains
¶ Get the number of chains in the Topology
-
property
n_residues
¶ Get the number of residues in the Topology.
-
residue
(self, index)¶ Get a specific residue by index. These indices start from zero.
- Parameters
- indexint
The index of the residue to select.
- Returns
- residueResidue
The index-th residue in the topology.
-
property
residues
¶ Iterator over all Residues in the Topology.
- Returns
- residueitergenerator
Iterator over all Residues in the Topology.
-
select
(self, selection_string)¶ Execute a selection against the topology
- Parameters
- selection_stringstr
An expression in the MDTraj atom selection DSL
- Returns
- indicesnp.ndarray, dtype=int, ndim=1
Array of the indices of the atoms matching the selection expression.
See also
select_expression
,mdtraj.core.selection.parse_selection
Examples
>>> topology.select('name O and water') array([1, 3, 5, 10, ...])
-
select_atom_indices
(self, selection='minimal')¶ Get the indices of biologically-relevant groups by name.
- Parameters
- selection{‘all’, ‘alpha’, ‘minimal’, ‘heavy’, ‘water’}
What types of atoms to select.
all
All atoms
alpha
Protein residue alpha carbons
minimal
Keep the atoms in protein residues with names in {CA, CB, C, N, O}
heavy
All non-hydrogen protein atoms.
water
Water oxygen atoms
- Returns
- indicesnp.ndarray (N,)
An array of the indices of the selected atoms.
-
select_expression
(self, selection_string)¶ Translate a atom selection expression into a pure python expression.
- Parameters
- selection_stringstr
An expression in the MDTraj atom selection DSL
- Returns
- python_stringstr
A string containing a pure python expression, equivalent to the selection expression.
Examples
>>> topology.select_expression('name O and water') "[atom.index for atom in topology.atoms if ((atom.name == 'O') and atom.residue.is_water)]")
-
select_pairs
(self, selection1=None, selection2=None)¶ Generate unique pairs of atom indices.
If a selecton is a string, it will be resolved using the atom selection DSL, otherwise it is expected to be an array of atom indices.
- Parameters
- selection1str or array-like, shape=(n_indices, ), dtype=int
A selection for select() or an array of atom indices.
- selection2str or array-like, shape=(n_indices, ), dtype=int
A selection for select() or an array of atom indices.
- Returns
- pairsarray-like, shape=(n_pairs, 2), dtype=int
Each row gives the indices of two atoms.
-
subset
(self, atom_indices)¶ Create a new Topology from a subset of the atoms in an existing topology.
- Parameters
- atom_indicesarray_like
A list of the indices corresponding to the atoms in that you’d like to retain.
Notes
The existing topology will not be altered.
-
to_bondgraph
(self)¶ Create a NetworkX graph from the atoms and bonds in this topology
- Returns
- gnx.Graph
A graph whose nodes are the Atoms in this topology, and whose edges are the bonds
Notes
This method requires the NetworkX python package.
-
to_dataframe
(self)¶ Convert this topology into a pandas dataframe
- Returns
- atomspandas.DataFrame
The atoms in the topology, represented as a data frame.
- bondsnp.ndarray, shape=(n_bonds, 4), dtype=float, Optional
The bonds in this topology, represented as an n_bonds x 4 array indicating the two atom indices of the bond, the bond type, and bond order, cast to floats as the type is mapped from the classes to fractional. The atom indices and order are integers cast to float
-
to_fasta
(self, chain=None)¶ Convert this topology into FASTA string
- Parameters
- chainInteger, optional, default=None
If specified, will return the FASTA string for this chain in the Topology.
- Returns
- fastaString or list of Strings
A FASTA string for each chain specified.
-
to_openmm
(self, traj=None)¶ Convert this topology into OpenMM topology
- Parameters
- trajMDTraj.Trajectory, optional, default=None
If specified, use the first frame from this trajectory to set the unitcell information in the openmm topology.
- Returns
- topologysimtk.openmm.app.Topology
This topology, as an OpenMM topology