openff.toolkit.topology.Topology

class openff.toolkit.topology.Topology(other=None)[source]

A Topology is a chemical representation of a system containing one or more molecules appearing in a specified order.

Warning

This API is experimental and subject to change.

Examples

Import some utilities

>>> from openmm import app
>>> from openff.toolkit.tests.utils import get_data_file_path, get_packmol_pdb_file_path
>>> pdb_filepath = get_packmol_pdb_file_path('cyclohexane_ethanol_0.4_0.6')
>>> monomer_names = ('cyclohexane', 'ethanol')

Create a Topology object from a PDB file and sdf files defining the molecular contents

>>> from openff.toolkit import Molecule, Topology
>>> pdbfile = app.PDBFile(pdb_filepath)
>>> sdf_filepaths = [get_data_file_path(f'systems/monomers/{name}.sdf') for name in monomer_names]
>>> unique_molecules = [Molecule.from_file(sdf_filepath) for sdf_filepath in sdf_filepaths]
>>> topology = Topology.from_openmm(pdbfile.topology, unique_molecules=unique_molecules)

Create a Topology object from a PDB file and IUPAC names of the molecular contents

>>> pdbfile = app.PDBFile(pdb_filepath)
>>> unique_molecules = [Molecule.from_iupac(name) for name in monomer_names]
>>> topology = Topology.from_openmm(pdbfile.topology, unique_molecules=unique_molecules)

Create an empty Topology object and add a few copies of a single benzene molecule

>>> topology = Topology()
>>> molecule = Molecule.from_iupac('benzene')
>>> molecule_topology_indices = [topology.add_molecule(molecule) for index in range(10)]
__init__(other=None)[source]

Create a new Topology.

Parameters

other (optional, default=None) – If specified, attempt to construct a copy of the Topology from the specified object. This might be a Topology object, or a file that can be used to construct a Topology object or serialized Topology object.

Methods

__init__([other])

Create a new Topology.

add_constraint(iatom, jatom[, distance])

Mark a pair of atoms as constrained.

add_molecule(molecule)

assert_bonded(atom1, atom2)

Raise an exception if the specified atoms are not bonded in the topology.

atom(atom_topology_index)

Get the TopologyAtom at a given Topology atom index.

atom_index(atom)

Returns the index of a given atom in this topology

bond(bond_topology_index)

Get the TopologyBond at a given Topology bond index.

chemical_environment_matches(query[, ...])

Retrieve all matches for a given chemical environment query.

copy_initializer(other)

from_bson(serialized)

Instantiate an object from a BSON serialized representation.

from_dict(topology_dict)

Create a new Topology from a dictionary representation

from_json(serialized)

Instantiate an object from a JSON serialized representation.

from_mdtraj(mdtraj_topology[, unique_molecules])

Construct an OpenFF Topology object from an MDTraj Topology object.

from_messagepack(serialized)

Instantiate an object from a MessagePack serialized representation.

from_molecules(molecules)

Create a new Topology object containing one copy of each of the specified molecule(s).

from_openmm(openmm_topology[, unique_molecules])

Construct an OpenFF Topology object from an OpenMM Topology object.

from_pickle(serialized)

Instantiate an object from a pickle serialized representation.

from_toml(serialized)

Instantiate an object from a TOML serialized representation.

from_xml(serialized)

Instantiate an object from an XML serialized representation.

from_yaml(serialized)

Instantiate from a YAML serialized representation.

get_bond_between(i, j)

Returns the bond between two atoms

hierarchy_iterator(iter_name)

Get a HierarchyElement iterator from all of the molecules in this topology that provide the appropriately named iterator.

is_bonded(i, j)

Returns True if the two atoms are bonded

is_constrained(iatom, jatom)

Check if a pair of atoms are marked as constrained.

molecule(index)

Returns the molecule with a given index in this Topology.

molecule_atom_start_index(molecule)

Returns the index of a molecule's first atom in this topology

molecule_index(molecule)

Returns the index of a given molecule in this topology

nth_degree_neighbors(n_degrees)

Return canonicalized pairs of atoms whose shortest separation is exactly n bonds.

particle_index(particle)

DEPRECATED: Use Topology.atom_index instead.

to_bson()

Return a BSON serialized representation.

to_dict()

to_file(filename, positions[, file_format, ...])

Save coordinates and topology to a PDB file.

to_json([indent])

Return a JSON serialized representation.

to_messagepack()

Return a MessagePack representation.

to_openmm([ensure_unique_atom_names])

Create an OpenMM Topology object.

to_pickle()

Return a pickle serialized representation.

to_toml()

Return a TOML serialized representation.

to_xml([indent])

Return an XML representation.

to_yaml()

Return a YAML serialized representation.

Attributes

amber_impropers

Iterate over improper torsions in the molecule, but only those with trivalent centers, reporting the central atom first in each improper.

angles

iterator over the angles in this Topology.

aromaticity_model

Get the aromaticity model applied to all molecules in the topology.

atoms

Returns an iterator over the atoms in this Topology.

bonds

Returns an iterator over the bonds in this Topology

box_vectors

Return the box vectors of the topology, if specified

constrained_atom_pairs

Returns the constrained atom pairs of the Topology

identical_molecule_groups

Returns groups of chemically identical molecules, identified by index and atom map.

impropers

iterator over the possible improper torsions in this Topology.

is_periodic

Return whether or not this Topology is intended to be described with periodic boundary conditions.

molecules

Returns an iterator over all the Molecules in this Topology

n_angles

number of angles in this Topology.

n_atoms

Returns the number of atoms in in this Topology.

n_bonds

Returns the number of Bonds in in this Topology.

n_impropers

number of possible improper torsions in this Topology.

n_molecules

Returns the number of molecules in this Topology

n_particles

Use Topology.n_atoms instead.

n_propers

number of proper torsions in this Topology.

n_reference_molecules

Use Topology.n_molecules instead.

n_topology_atoms

Use Topology.n_atoms instead.

n_topology_bonds

Use Topology.n_bonds instead.

n_topology_molecules

Use Topology.n_molecules instead.

n_topology_particles

Use Topology.n_particles instead.

particles

Use Topology.molecules instead.

propers

iterator over the proper torsions in this Topology.

reference_molecules

Get a list of reference molecules in this Topology.

smirnoff_impropers

Iterate over improper torsions in the molecule, but only those with trivalent centers, reporting the central atom second in each improper.

topology_atoms

Use Topology.atoms instead.

topology_bonds

Use Topology.bonds instead.

topology_molecules

Use Topology.molecules instead.

topology_particles

Use Topology.particles instead.

property reference_molecules: List[Molecule]

Get a list of reference molecules in this Topology.

Returns

iterable of openff.toolkit.topology.Molecule

classmethod from_molecules(molecules: Union[Molecule, List[Molecule]])[source]

Create a new Topology object containing one copy of each of the specified molecule(s).

Parameters

molecules (Molecule or iterable of Molecules) – One or more molecules to be added to the Topology

Returns

topology (Topology) – The Topology created from the specified molecule(s)

assert_bonded(atom1, atom2)[source]

Raise an exception if the specified atoms are not bonded in the topology.

Parameters
  • atom1 (Atom or int) – The atoms or atom topology indices to check to ensure they are bonded

  • atom2 (Atom or int) – The atoms or atom topology indices to check to ensure they are bonded

property aromaticity_model

Get the aromaticity model applied to all molecules in the topology.

Returns

aromaticity_model (str) – Aromaticity model in use.

property box_vectors

Return the box vectors of the topology, if specified

Returns

box_vectors (unit-wrapped numpy array of shape (3, 3)) – The unit-wrapped box vectors of this topology

property is_periodic

Return whether or not this Topology is intended to be described with periodic boundary conditions.

property constrained_atom_pairs: Dict[Tuple[int], Union[openff.units.units.Quantity, bool]]

Returns the constrained atom pairs of the Topology

Returns

constrained_atom_pairs (dict of Tuple[int]: Union[unit.Quantity, bool]) – dictionary of the form {(atom1_topology_index, atom2_topology_index): distance}

property n_molecules

Returns the number of molecules in this Topology

Returns

n_molecules (Iterable of Molecule)

property molecules: Generator[Union[Molecule, openff.toolkit.topology._mm_molecule._SimpleMolecule], None, None]

Returns an iterator over all the Molecules in this Topology

Returns

molecules (Iterable of Molecule)

molecule(index)[source]

Returns the molecule with a given index in this Topology.

Returns

molecule (openff.toolkit.topology.Molecule)

property n_atoms: int

Returns the number of atoms in in this Topology.

Returns

n_atoms (int)

property atoms: Generator[Atom, None, None]

Returns an iterator over the atoms in this Topology. These will be in ascending order of topology index.

Returns

atoms (Generator of TopologyAtom)

atom_index(atom)[source]

Returns the index of a given atom in this topology

Parameters

atom (Atom) –

Returns

index (int) – The index of the given atom in this topology

:raises AtomNotInTopologyError : If the given atom is not in this topology:

molecule_index(molecule)[source]

Returns the index of a given molecule in this topology

Parameters

molecule (FrozenMolecule) –

Returns

index (int) – The index of the given molecule in this topology

:raises MoleculeNotInTopologyError : If the given atom is not in this topology:

molecule_atom_start_index(molecule)[source]

Returns the index of a molecule’s first atom in this topology

Parameters

molecule (FrozenMolecule) –

Returns

index (int)

property n_bonds

Returns the number of Bonds in in this Topology.

Returns

n_bonds (int)

property bonds

Returns an iterator over the bonds in this Topology

Returns

bonds (Iterable of Bond)

property n_angles

number of angles in this Topology.

Type

int

property angles

iterator over the angles in this Topology.

Type

Iterable of Tuple[Atom]

property n_propers

number of proper torsions in this Topology.

Type

int

property propers

iterator over the proper torsions in this Topology.

Type

Iterable of Tuple[TopologyAtom]

property n_impropers

number of possible improper torsions in this Topology.

Type

int

property impropers

iterator over the possible improper torsions in this Topology.

Type

Iterable of Tuple[TopologyAtom]

property smirnoff_impropers

Iterate over improper torsions in the molecule, but only those with trivalent centers, reporting the central atom second in each improper.

Note that it’s possible that a trivalent center will not have an improper assigned. This will depend on the force field that is used.

Also note that this will return 6 possible atom orderings around each improper center. In current SMIRNOFF parameterization, three of these six orderings will be used for the actual assignment of the improper term and measurement of the angles. These three orderings capture the three unique angles that could be calculated around the improper center, therefore the sum of these three terms will always return a consistent energy.

For more details on the use of three-fold (‘trefoil’) impropers, see https://openforcefield.github.io/standards/standards/smirnoff/#impropertorsions

Returns

impropers (set of tuple) – An iterator of tuples, each containing the indices of atoms making up a possible improper torsion. The central atom is listed second in each tuple.

property amber_impropers

Iterate over improper torsions in the molecule, but only those with trivalent centers, reporting the central atom first in each improper.

Note that it’s possible that a trivalent center will not have an improper assigned. This will depend on the force field that is used.

Also note that this will return 6 possible atom orderings around each improper center. In current AMBER parameterization, one of these six orderings will be used for the actual assignment of the improper term and measurement of the angle. This method does not encode the logic to determine which of the six orderings AMBER would use.

Returns

impropers (set of tuple) – An iterator of tuples, each containing the indices of atoms making up a possible improper torsion. The central atom is listed first in each tuple.

nth_degree_neighbors(n_degrees: int)[source]

Return canonicalized pairs of atoms whose shortest separation is exactly n bonds. Only pairs with increasing atom indices are returned.

Parameters

n (int) – The number of bonds separating atoms in each pair

Returns

neighbors (iterator of tuple of Atom) – Tuples (len 2) of atom that are separated by n bonds.

Notes

The criteria used here relies on minimum distances; when there are multiple valid paths between atoms, such as atoms in rings, the shortest path is considered. For example, two atoms in “meta” positions with respect to each other in a benzene are separated by two paths, one length 2 bonds and the other length 4 bonds. This function would consider them to be 2 apart and would not include them if n=4 was passed.

chemical_environment_matches(query, aromaticity_model='MDL', unique=False, toolkit_registry=GLOBAL_TOOLKIT_REGISTRY)[source]

Retrieve all matches for a given chemical environment query.

TODO:

  • Do we want to generalize this to other kinds of queries too, like mdtraj DSL, pymol selections, atom index slices, etc? We could just call it topology.matches(query)

Parameters
  • query (str or ChemicalEnvironment) – SMARTS string (with one or more tagged atoms) or ChemicalEnvironment query Query will internally be resolved to SMARTS using query.as_smarts() if it has an .as_smarts method.

  • aromaticity_model (str) – Override the default aromaticity model for this topology and use the specified aromaticity model instead. Allowed values: [‘MDL’]

Returns

matches (list of Topology._ChemicalEnvironmentMatch) – A list of tuples, containing the topology indices of the matching atoms.

property identical_molecule_groups

Returns groups of chemically identical molecules, identified by index and atom map.

Returns

  • identical_molecule_groups ({int:[[int: {int: int}]]}) – A dict of the form {unique_mol_idx : [[topology_mol_idx, atom_map],…]. A dict where each key is the topology molecule index of a unique chemical species, and each value is a list describing all of the instances of that chemical species in the topology. Each instance is a two-membered list where the first element is the topology molecule index, and the second element is a dict describing the atom map from the unique molecule to the instance of it in the topology.

  • >>> from openff.toolkit import Molecule, Topology

  • >>> # Create a water ordered as OHH

  • >>> water1 = Molecule()

  • >>> water1.add_atom(8, 0, False)

  • >>> water1.add_atom(1, 0, False)

  • >>> water1.add_atom(1, 0, False)

  • >>> water1.add_bond(0, 1, 1, False)

  • >>> water1.add_bond(0, 2, 1, False)

  • >>> # Create a different water ordered as HOH

  • >>> water2 = Molecule()

  • >>> water2.add_atom(1, 0, False)

  • >>> water2.add_atom(8, 0, False)

  • >>> water2.add_atom(1, 0, False)

  • >>> water2.add_bond(0, 1, 1, False)

  • >>> water2.add_bond(1, 2, 1, False)

  • >>> top = Topology.from_molecules([water1, water2])

  • >>> top.identical_molecule_groups

  • {0 ([[0, {0: 0, 1: 1, 2: 2}], [1, {0: 1, 1: 0, 2: 2}]]})

classmethod from_dict(topology_dict)[source]

Create a new Topology from a dictionary representation

Parameters

topology_dict (OrderedDict) – A dictionary representation of the topology.

Returns

topology (Topology) – A Topology created from the dictionary representation

classmethod from_openmm(openmm_topology, unique_molecules=None)[source]

Construct an OpenFF Topology object from an OpenMM Topology object.

This method guarantees that the order of atoms in the input OpenMM Topology will be the same as the ordering of atoms in the output OpenFF Topology. However it does not guarantee the order of the bonds will be the same.

Parameters
  • openmm_topology (openmm.app.Topology) – An OpenMM Topology object

  • unique_molecules (iterable of objects that can be used to construct unique Molecule objects) – All unique molecules must be provided, in any order, though multiple copies of each molecule are allowed. The atomic elements and bond connectivity will be used to match the reference molecules to molecule graphs appearing in the OpenMM Topology. If bond orders are present in the OpenMM Topology, these will be used in matching as well.

Returns

topology (openff.toolkit.topology.Topology) – An OpenFF Topology object

to_openmm(ensure_unique_atom_names=True)[source]

Create an OpenMM Topology object.

The atom metadata fields residue_name, residue_number, and chain_id are used to group atoms into OpenMM residues and chains.

Contiguously-indexed atoms with the same residue_name, residue_number, and chain_id will be put into the same OpenMM residue.

Contiguously-indexed residues with with the same chain_id will be put into the same OpenMM chain.

This method will never make an OpenMM chain or residue larger than the OpenFF Molecule that it came from. In other words, no chain or residue will span two OpenFF Molecules.

This method will not populate the OpenMM Topology with virtual sites.

Parameters

ensure_unique_atom_names (bool, optional. Default=True) – Whether to check that the molecules in each molecule have unique atom names, and regenerate them if not. Note that this looks only at molecules, and does not guarantee uniqueness in the entire Topology.

Returns

openmm_topology (openmm.app.Topology) – An OpenMM Topology object

to_file(filename: str, positions, file_format: str = 'PDB', keepIds=False)[source]

Save coordinates and topology to a PDB file.

Reference: https://github.com/openforcefield/openff-toolkit/issues/502

Notes:

  1. Atom numbering may not remain same, for example if the atoms in water are numbered as 1001, 1002, 1003, they would change to 1, 2, 3. This doesn’t affect the topology or coordinates or atom-ordering in any way.

  2. Same issue with the amino acid names in the pdb file, they are not returned.

Parameters
  • filename (str) – name of the pdb file to write to

  • positions (n_atoms x 3 numpy array or openmm.unit.Quantity-wrapped n_atoms x 3 iterable) –

    Can be a
    • openmm.unit.Quantity’ object which has atomic positions as a list of unit-tagged `Vec3 objects

    • openff.units.unit.Quantity object which wraps a numpy.ndarray with units

    • (unitless) 2D numpy.ndarray, in which it is assumed that the positions are in units of Angstroms.

    For all data types, must have shape (n_atoms, 3) where n_atoms is the number of atoms in this topology.

  • file_format (str) – Output file format. Case insensitive. Currently only supported value is “pdb”.

classmethod from_mdtraj(mdtraj_topology, unique_molecules=None)[source]

Construct an OpenFF Topology object from an MDTraj Topology object.

Parameters
  • mdtraj_topology (Topology) – An MDTraj Topology object

  • unique_molecules (iterable of objects that can be used to construct unique Molecule objects) – All unique molecules must be provided, in any order, though multiple copies of each molecule are allowed. The atomic elements and bond connectivity will be used to match the reference molecules to molecule graphs appearing in the MDTraj Topology. If bond orders are present in the MDTraj Topology, these will be used in matching as well.

Returns

topology (openff.toolkit.topology.Topology) – An OpenFF Topology object

get_bond_between(i, j)[source]

Returns the bond between two atoms

Parameters
  • i (int or TopologyAtom) – Atoms or atom indices to check

  • j (int or TopologyAtom) – Atoms or atom indices to check

Returns

bond (TopologyBond) – The bond between i and j.

is_bonded(i, j)[source]

Returns True if the two atoms are bonded

Parameters
  • i (int or TopologyAtom) – Atoms or atom indices to check

  • j (int or TopologyAtom) – Atoms or atom indices to check

Returns

is_bonded (bool) – True if atoms are bonded, False otherwise.

atom(atom_topology_index)[source]

Get the TopologyAtom at a given Topology atom index.

Parameters

atom_topology_index (int) – The index of the TopologyAtom in this Topology

Returns

An openff.toolkit.topology.TopologyAtom

bond(bond_topology_index)[source]

Get the TopologyBond at a given Topology bond index.

Parameters

bond_topology_index (int) – The index of the TopologyBond in this Topology

Returns

An openff.toolkit.topology.TopologyBond

add_constraint(iatom, jatom, distance=True)[source]

Mark a pair of atoms as constrained.

Constraints between atoms that are not bonded (e.g., rigid waters) are permissible.

Parameters
  • iatom (Atom) – Atoms to mark as constrained These atoms may be bonded or not in the Topology

  • jatom (Atom) – Atoms to mark as constrained These atoms may be bonded or not in the Topology

  • distance (unit-wrapped float, optional, default=True) – Constraint distance True if distance has yet to be determined False if constraint is to be removed

is_constrained(iatom, jatom)[source]

Check if a pair of atoms are marked as constrained.

Parameters
  • iatom (int) – Indices of atoms to mark as constrained.

  • jatom (int) – Indices of atoms to mark as constrained.

Returns

distance (unit-wrapped float or bool) – True if constrained but constraints have not yet been applied Distance if constraint has already been added to Topology

hierarchy_iterator(iter_name)[source]

Get a HierarchyElement iterator from all of the molecules in this topology that provide the appropriately named iterator. This iterator will yield HierarchyElements sorted first by the order that molecules are listed in the Topology, and second by the specific sorting of HierarchyElements defined in each molecule.

Parameters

iter_name (string) – The iterator name associated with the HierarchyScheme to retrieve (for example ‘residues’ or ‘chains’)

Returns

iterator of HierarchyElement

property n_topology_atoms: int

Use Topology.n_atoms instead.

Type

DEPRECATED

property topology_atoms

Use Topology.atoms instead.

Type

DEPRECATED

property n_topology_bonds: int

Use Topology.n_bonds instead.

Type

DEPRECATED

property topology_bonds

Use Topology.bonds instead.

Type

DEPRECATED

property n_topology_particles: int

Use Topology.n_particles instead.

Type

DEPRECATED

property topology_particles

Use Topology.particles instead.

Type

DEPRECATED

property n_reference_molecules: int

Use Topology.n_molecules instead.

Type

DEPRECATED

property n_topology_molecules: int

Use Topology.n_molecules instead.

Type

DEPRECATED

property topology_molecules

Use Topology.molecules instead.

Type

DEPRECATED

classmethod from_bson(serialized)

Instantiate an object from a BSON serialized representation.

Specification: http://bsonspec.org/

Parameters

serialized (bytes) – A BSON serialized representation of the object

Returns

instance (cls) – An instantiated object

classmethod from_json(serialized: str)

Instantiate an object from a JSON serialized representation.

Specification: https://www.json.org/

Parameters

serialized (str) – A JSON serialized representation of the object

Returns

instance (cls) – An instantiated object

classmethod from_messagepack(serialized)

Instantiate an object from a MessagePack serialized representation.

Specification: https://msgpack.org/index.html

Parameters

serialized (bytes) – A MessagePack-encoded bytes serialized representation

Returns

instance (cls) – Instantiated object.

classmethod from_pickle(serialized)

Instantiate an object from a pickle serialized representation.

Warning

This is not recommended for safe, stable storage since the pickle specification may change between Python versions.

Parameters

serialized (str) – A pickled representation of the object

Returns

instance (cls) – An instantiated object

classmethod from_toml(serialized)

Instantiate an object from a TOML serialized representation.

Specification: https://github.com/toml-lang/toml

Parameters

serlialized (str) – A TOML serialized representation of the object

Returns

instance (cls) – An instantiated object

classmethod from_xml(serialized)

Instantiate an object from an XML serialized representation.

Specification: https://www.w3.org/XML/

Parameters

serialized (bytes) – An XML serialized representation

Returns

instance (cls) – Instantiated object.

classmethod from_yaml(serialized)

Instantiate from a YAML serialized representation.

Specification: http://yaml.org/

Parameters

serialized (str) – A YAML serialized representation of the object

Returns

instance (cls) – Instantiated object

property n_particles: int

Use Topology.n_atoms instead.

Type

DEPRECATED

to_bson()

Return a BSON serialized representation.

Specification: http://bsonspec.org/

Returns

serialized (bytes) – A BSON serialized representation of the objecft

to_json(indent=None) str

Return a JSON serialized representation.

Specification: https://www.json.org/

Parameters

indent (int, optional, default=None) – If not None, will pretty-print with specified number of spaces for indentation

Returns

serialized (str) – A JSON serialized representation of the object

to_messagepack()

Return a MessagePack representation.

Specification: https://msgpack.org/index.html

Returns

serialized (bytes) – A MessagePack-encoded bytes serialized representation of the object

to_pickle()

Return a pickle serialized representation.

Warning

This is not recommended for safe, stable storage since the pickle specification may change between Python versions.

Returns

serialized (str) – A pickled representation of the object

to_toml()

Return a TOML serialized representation.

Specification: https://github.com/toml-lang/toml

Returns

serialized (str) – A TOML serialized representation of the object

to_xml(indent=2)

Return an XML representation.

Specification: https://www.w3.org/XML/

Parameters

indent (int, optional, default=2) – If not None, will pretty-print with specified number of spaces for indentation

Returns

serialized (bytes) – A MessagePack-encoded bytes serialized representation.

to_yaml()

Return a YAML serialized representation.

Specification: http://yaml.org/

Returns

serialized (str) – A YAML serialized representation of the object

property particles

Use Topology.molecules instead.

Type

DEPRECATED

particle_index(particle) int[source]

DEPRECATED: Use Topology.atom_index instead.