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
chains | Iterator over all Chains in the Topology. |
residues | Iterator over all Residues in the Topology. |
atoms | Iterator over all Atoms in the Topology. |
Create a new Topology object
Methods
__init__() | Create a new Topology object |
add_atom(name, element, residue[, serial]) | Create a new Atom and add it to the Topology. |
add_bond(atom1, atom2) | Create a new bond and add it to the Topology. |
add_chain() | Create a new Chain and add it to the Topology. |
add_residue(name, chain[, resSeq]) | Create a new Residue and add it to the Topology. |
atom(index) | Get a specific atom by index. |
atoms_by_name(name) | Iterator over all Atoms in the Topology with a specified name |
chain(index) | Get a specific chain by index. |
copy() | Return a copy of the topology |
create_disulfide_bonds(positions) | Identify disulfide bonds based on proximity and add them to the Topology. |
create_standard_bonds() | Create bonds based on the atom and residue names for all standard residue types. |
from_dataframe(atoms[, bonds]) | Create a mdtraj topology from a pandas data frame |
from_openmm(value) | Create a mdtraj topology from an OpenMM topology |
join(other) | Join two topologies together |
residue(index) | Get a specific residue by index. |
select(selection_string) | Execute a selection against the topology |
select_atom_indices([selection]) | Get the indices of biologically-relevant groups by name. |
select_expression(selection_string) | Translate a atom selection expression into a pure python expression. |
select_pairs([selection1, selection2]) | Generate unique pairs of atom indices. |
subset(atom_indices) | Create a new Topology from a subset of the atoms in an existing topology. |
to_bondgraph() | Create a NetworkX graph from the atoms and bonds in this topology |
to_dataframe() | Convert this topology into a pandas dataframe |
to_fasta([chain]) | Convert this topology into FASTA string |
to_openmm([traj]) | Convert this topology into OpenMM topology |
Attributes
atoms | Iterator over all Atoms in the Topology. |
bonds | Iterator over all bonds (each represented as a tuple of two Atoms) in the Topology. |
chains | Iterator over all Chains in the Topology. |
n_atoms | Get the number of atoms in the Topology |
n_bonds | Get the number of bonds in the Topology |
n_chains | Get the number of chains in the Topology |
n_residues | Get the number of residues in the Topology. |
residues | Iterator over all Residues in the Topology. |
Return a copy of the topology
Returns: | out : Topology
|
---|
Join two topologies together
Parameters: | other : Topology
|
---|---|
Returns: | out : Topology
|
Convert this topology into FASTA string
Parameters: | chain : Integer, optional, default=None
|
---|---|
Returns: | fasta : String or list of Strings
|
Convert this topology into OpenMM topology
Parameters: | traj : MDTraj.Trajectory, optional, default=None
|
---|---|
Returns: | topology : simtk.openmm.app.Topology
|
Create a mdtraj topology from an OpenMM topology
Parameters: | value : simtk.openmm.app.Topology
|
---|
Convert this topology into a pandas dataframe
Returns: | atoms : pandas.DataFrame
bonds : np.ndarray
|
---|
Create a mdtraj topology from a pandas data frame
Parameters: | atoms : pandas.DataFrame
bonds : np.ndarray, shape=(n_bonds, 2), dtype=int, optional
|
---|
See also
Create a NetworkX graph from the atoms and bonds in this topology
Returns: | g : nx.Graph
|
---|
Notes
This method requires the NetworkX python package.
Create a new Chain and add it to the Topology.
Returns: | chain : mdtraj.topology.Chain
|
---|
Create a new Residue and add it to the Topology.
Parameters: | name : str
chain : mdtraj.topology.Chain
resSeq : int, optional
|
---|---|
Returns: | residue : mdtraj.topology.Residue
|
Create a new Atom and add it to the Topology.
Parameters: | name : str
element : mdtraj.element.Element
residue : mdtraj.topology.Residue
serial : int
|
---|---|
Returns: | atom : mdtraj.topology.Atom
|
Create a new bond and add it to the Topology.
Parameters: | atom1 : mdtraj.topology.Atom
atom2 : mdtraj.topology.Atom
|
---|
Get a specific chain by index. These indices start from zero.
Parameters: | index : int
|
---|---|
Returns: | chain : Chain
|
Iterator over all Chains in the Topology.
Returns: | chainiter : listiterator
|
---|
Get the number of chains in the Topology
Get a specific residue by index. These indices start from zero.
Parameters: | index : int
|
---|---|
Returns: | residue : Residue
|
Iterator over all Residues in the Topology.
Returns: | residueiter : generator
|
---|
Get the number of residues in the Topology.
Get a specific atom by index. These indices start from zero.
Parameters: | index : int
|
---|---|
Returns: | atom : Atom
|
Iterator over all Atoms in the Topology.
Returns: | atomiter : generator
|
---|
Iterator over all Atoms in the Topology with a specified name
Parameters: | name : str
|
---|---|
Returns: | atomiter : generator |
Examples
>>> for atom in topology.atoms_by_name('CA'):
... print(atom)
Get the number of atoms in the Topology
Iterator over all bonds (each represented as a tuple of two Atoms) in the Topology.
Returns: | atomiter : generator
|
---|
Get the number of bonds in the Topology
Create bonds based on the atom and residue names for all standard residue types.
Identify disulfide bonds based on proximity and add them to the Topology.
Parameters: | positions : list
|
---|
Create a new Topology from a subset of the atoms in an existing topology.
Parameters: | atom_indices : array_like
|
---|
Notes
The existing topology will not be altered.
Translate a atom selection expression into a pure python expression.
Parameters: | selection_string : str
|
---|---|
Returns: | python_string : str
|
Examples
>>> topology.select_expression('name O and water')
"[atom.index for atom in topology.atoms if ((atom.name == 'O') and atom.residue.is_water)]")
Execute a selection against the topology
Parameters: | selection_string : str
|
---|---|
Returns: | indices : np.ndarray, dtype=int, ndim=1
|
See also
select_expression, mdtraj.core.selection.parse_selection
Examples
>>> topology.select('name O and water')
array([1, 3, 5, 10, ...])
Get the indices of biologically-relevant groups by name.
Parameters: | selection : {‘all’, ‘alpha’, ‘minimal’, ‘heavy’, ‘water’}
|
---|---|
Returns: | indices : np.ndarray (N,)
|
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: | selection1 : str or array-like, shape=(n_indices, ), dtype=int
selection2 : str or array-like, shape=(n_indices, ), dtype=int
|
---|---|
Returns: | pairs : array-like, shape=(n_pairs, 2), dtype=int
|