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

Add a copy of the molecule to the topology

assert_bonded(atom1, atom2)

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

atom(atom_topology_index)

Get the Atom 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 Bond at a given Topology bond index.

chemical_environment_matches(query[, ...])

Retrieve all matches for a given chemical environment query.

clear_positions()

Clear the positions of this topology by removing all conformers from its consituent molecules.

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[, ...])

Construct an OpenFF Topology from an MDTraj Topology

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[, ...])

Construct an OpenFF Topology object from an OpenMM Topology object.

from_pdb(file_path[, unique_molecules, ...])

Loads supported or user-specified molecules from a PDB file.

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

get_positions()

Copy the positions of the topology into a new array.

hierarchy_iterator(iter_name)

Iterate over all molecules with the given hierarchy scheme.

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.

set_positions(array)

Set the positions in a topology by copying from a single n×3 array.

to_bson()

Return a BSON serialized representation.

to_dict()

to_file(file[, 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.

visualize([ensure_correct_connectivity])

Visualize with NGLView.

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

The number of possible improper torsions in this Topology.

n_molecules

Returns the number of molecules in this Topology

n_propers

The number of proper torsions in this Topology.

n_unique_molecules

Returns the number of unique molecules in this Topology

propers

iterator over the proper torsions 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.

unique_molecules

Get a list of chemically unique molecules in this Topology.

property unique_molecules: Iterator[Molecule | _SimpleMolecule]

Get a list of chemically unique molecules in this Topology.

Molecules are considered unique if they fail an isomorphism check with default values (see Molecule.is_isomorphic_with()). The order of molecules returned by this property is arbitrary.

property n_unique_molecules: int

Returns the number of unique molecules in this Topology

classmethod from_molecules(molecules: Molecule | FrozenMolecule | _SimpleMolecule | list[Union[Molecule, FrozenMolecule, openff.toolkit.topology._mm_molecule._SimpleMolecule]]) Topology[source]

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

Parameters:

molecules – One or more molecules to be added to the Topology

Returns:

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

assert_bonded(atom1: int | Atom, atom2: int | Atom)[source]

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

Parameters:
  • atom1 – The atoms or atom topology indices to check to ensure they are bonded

  • atom2 – The atoms or atom topology indices to check to ensure they are bonded

property aromaticity_model: str

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

Returns:

aromaticity_model – Aromaticity model in use.

property box_vectors

Return the box vectors of the topology, if specified

Returns:

box_vectors – The unit-wrapped box vectors of this topology

property is_periodic: bool

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

property constrained_atom_pairs: dict[tuple[int, int], Union[Quantity, bool]]

Returns the constrained atom pairs of the Topology

Returns:

constrained_atom_pairs – dictionary of the form {(atom1_topology_index, atom2_topology_index): distance}

property n_molecules: int

Returns the number of molecules in this Topology

property molecules: Generator[Molecule | FrozenMolecule | _SimpleMolecule, None, None]

Returns an iterator over all the Molecules in this Topology

Returns:

molecules

molecule(index: int) Molecule | _SimpleMolecule[source]

Returns the molecule with a given index in this Topology.

Returns:

molecule

property n_atoms: int

Returns the number of atoms in in this Topology.

Returns:

n_atoms

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

atom_index(atom: Atom) int[source]

Returns the index of a given atom in this topology

Parameters:

atom

Returns:

index – The index of the given atom in this topology

Raises:

AtomNotInTopologyError – If the given atom is not in this topology

molecule_index(molecule: Molecule | FrozenMolecule | _SimpleMolecule) int[source]

Returns the index of a given molecule in this topology

Parameters:

molecule

Returns:

index – 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: Molecule) int[source]

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

Parameters:

molecule

Returns:

index

property n_bonds: int

Returns the number of Bonds in in this Topology.

Returns:

n_bonds

property bonds: Generator[Bond, None, None]

Returns an iterator over the bonds in this Topology

Returns:

bonds

property n_angles: int

number of angles in this Topology.

Type:

int

property angles: Generator[tuple['Atom', ...], None, None]

Iterator over the angles in this Topology. Returns a Generator of tuple[Atom].

property n_propers: int

The number of proper torsions in this Topology.

property propers: ], None, None]

iterator over the proper torsions in this Topology.

Type:

Iterable of tuple[Atom]

property n_impropers: int

The number of possible improper torsions in this Topology.

property impropers: Generator[tuple['Atom', ...], None, None]

iterator over the possible improper torsions in this Topology.

Type:

Generator of tuple[Atom]

property smirnoff_impropers: ], None, None]

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:

smirnoff_impropers – An iterator of tuples, each containing the Atom objects comprising up a possible improper torsion. The central atom is listed second in each tuple.

property amber_impropers: ], None, None]

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:

amber_impropers – An iterator of tuples, each containing the Atom objects comprising 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: str, aromaticity_model: str = 'MDL', unique: bool = False, toolkit_registry: ToolkitRegistry | ToolkitWrapper = 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 – SMARTS string (with one or more tagged atoms)

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

Returns:

matches – A list of tuples, containing the topology indices of the matching atoms.

property identical_molecule_groups: dict[int, list[tuple[int, dict[int, int]]]]

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

Returns:

identical_molecule_groups – A dict of the form {unique_mol_idx : [(topology_mol_idx, atom_map), ...]. Each key is the topology molecule index of the first instance of a unique chemical species. Iterating over the keys will yield all of the unique chemical species in the topology. Each value is a list of all instances of that chemical species in the topology. Each element of the list is a 2-tuple where the first element is the molecule index of the instance, and the second element maps the atom topology indices of the key molecule to the instance. The molecule instance corresponding to the key is included in the list, so the list is a complete list of all instances of that chemical species.

Examples

>>> from openff.toolkit import Molecule, Topology
>>> # Create a water ordered as OHH
>>> water1 = Molecule()
>>> water1.add_atom(8, 0, False)
0
>>> water1.add_atom(1, 0, False)
1
>>> water1.add_atom(1, 0, False)
2
>>> _ = 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)
0
>>> water2.add_atom(8, 0, False)
1
>>> water2.add_atom(1, 0, False)
2
>>> _ = 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: dict)[source]

Create a new Topology from a dictionary representation

Parameters:

topology_dict – A dictionary representation of the topology.

Returns:

topology – A Topology created from the dictionary representation

classmethod from_openmm(openmm_topology: openmm.app.Topology, unique_molecules: Iterable[FrozenMolecule] | None = None, positions: None | Quantity | OMMQuantity = None) Topology[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.

Hierarchy schemes are taken from the OpenMM topology, not from unique_molecules.

If any virtual sites are detected in the OpenMM topology, VirtualSitesUnsupportedError is raised because the Topology object model does not store virtual sites.

Parameters:
  • openmm_topology – The OpenMM Topology object to convert

  • unique_molecules – An iterable containing all the unique molecules in the topology. This is used to identify the molecules in the OpenMM topology and provide any missing chemical information. Each chemical species in the topology must be specified exactly once, though the topology may have any number of copies, including zero. The chemical elements of atoms and their bond connectivity will be used to match these reference molecules to the molecules appearing in the topology. If bond orders are specified in the topology, these will be used in matching as well.

  • positions – Positions for the atoms in the new topology.

Returns:

topology – An OpenFF Topology object, constructed from the molecules in unique_molecules, with the same atom order as the input topology.

Raises:
classmethod from_pdb(file_path: str | Path | TextIO, unique_molecules: Iterable[Molecule] | None = None, toolkit_registry=GLOBAL_TOOLKIT_REGISTRY, _custom_substructures: dict[str, list[str]] | None = None, _additional_substructures: Iterable[Molecule] | None = None)[source]

Loads supported or user-specified molecules from a PDB file.

Currently, canonical proteins, waters, and monoatomic ions are supported without CONECT records via residue and atom names, and molecules specified in the unique_molecules argument are supported when CONECT records are provided.

Warning

Molecules in the resulting Topology will adopt the geometric stereochemistry in the PDB, even if this conflicts with the stereochemistry specified in unique_molecules.

All molecules in the PDB file have the following requirements:

  • Polymer molecules must use the standard atom names described in the PDB Chemical Component Dictionary (PDB CCD).

  • There must be no missing atoms (all hydrogens must be explicit).

  • All particles must correspond to an atomic nucleus (particles in the PDB representing “virtual sites” or “extra points” are not allowed).

  • All bonds must be specified by either CONECT records, or for polymers and water by the PDB CCD via the residue and atom name.

  • CONECT records must correspond only to chemical bonds (CONECT records representing an angle constraints are not allowed).

  • CONECT records may be redundant with connectivity defined by residue templates.

Currently, the only supported polymers are proteins made of the 20 canonical amino acids. For details on the polymer loading used here, see Molecule.from_polymer_pdb().

Waters can be recognized in either of two ways:

  • By the residue name “HOH” and atom names “H1”, “H2”, and “O”.

  • By ATOM records which include element information and CONECT records.

Monoatomic ions supported by Sage are recognized (Na+, Li+, K+, Rb+, Cs+, F-, Cl-, Br-, and I-). To load other monoatomic ions, use the unique_molecules keyword argument.

The unique_molecules keyword argument can be used to load arbitrary molecules from the PDB file. These molecules match a group of atoms in the PDB file when their atomic elements and connectivity are identical; elements and CONECT records must therefore be explicitly specified in the PDB file. Information missing from the PDB format, such as bond orders and formal charges, is then taken from the matching unique molecule. Unique molecule matches will overwrite bond order and formal charge assignments from other sources. Stereochemistry is assigned based on the geometry in the PDB, even if this differs from the stereochemistry in the unique molecule.

A user-defined molecule in the PDB file must exactly match a unique molecule to successfully load it - substructures and superstructures will raise UnassignedChemistryInPDBError. Unique molecules need not be present in the PDB.

Metadata such as residues, chains, and atom names are recorded in the Atom.metadata attribute, which is a dictionary mapping from the strings "residue_name", "residue_number", "insertion_code", and "chain_id" to the appropriate value. The topology returned by this method can expose residue and chain iterators which can be accessed using Topology.hierarchy_iterator(), such as top.hierarchy_iterator("residues") and top.hierarchy_iterator ("chains").

A CRYST line in the PDB, if present, will be interpreted as periodic box vectors in Angstroms.

Parameters:
  • file_path – PDB information to be passed to OpenMM PDBFile object for loading

  • unique_molecules – OpenFF Molecule objects corresponding to the molecules in the input PDB. See above for details.

  • toolkit_registry – The ToolkitRegistry to use as the cheminformatics backend.

  • _custom_substructures (dict[str, list[str]], Default = {}) – Experimental and unstable. Dictionary where keys are the names of new substructures (cannot overlap with existing amino acid names) and the values are the new substructure entries that follow the same format as those used in the amino acid substructure library

  • _additional_substructures – Experimental and unstable. Molecule with atom.metadata[“substructure_atom”] = True or False for all atoms. Currently only stable for independent, standalone molecules not bonded to a larger protein/molecule. (For that use _custom_substructures)

Returns:

topology

Raises:

UnassignedChemistryInPDBError – If an atom or bond could not be assigned; the exception will provide a detailed diagnostic of what went wrong.

Examples

>>> from openff.toolkit import Topology
>>> from openff.toolkit.utils import get_data_file_path
>>> top = Topology.from_pdb(get_data_file_path("proteins/TwoMol_SER_CYS.pdb"))
>>> # The molecules in the loaded topology are full-fledged OpenFF Molecule objects
>>> for match in top.chemical_environment_matches('[O:1]=[C:2][N:3][H:4]'): print(match.topology_atom_indices)
(1, 0, 6, 13)
(9, 8, 17, 19)
(24, 23, 29, 36)
(32, 31, 40, 42)
>>> [*top.hierarchy_iterator("residues")]
[HierarchyElement ('A', '1', ' ', 'ACE') of iterator 'residues' containing 6 atom(s),
 HierarchyElement ('A', '2', ' ', 'SER') of iterator 'residues' containing 11 atom(s),
 HierarchyElement ('A', '3', ' ', 'NME') of iterator 'residues' containing 6 atom(s),
 HierarchyElement ('B', '1', ' ', 'ACE') of iterator 'residues' containing 6 atom(s),
 HierarchyElement ('B', '2', ' ', 'CYS') of iterator 'residues' containing 11 atom(s),
 HierarchyElement ('B', '3', ' ', 'NME') of iterator 'residues' containing 6 atom(s)]

Polymer systems can also be supported if _custom_substructures are given as a dict[str, list[str]], where the keys are unique atom names and the values are lists of substructure SMARTS. The substructure SMARTS must follow the same format as given in the residue substructure connectivity library: "<bond>[#<atomic number>D<degree>+<formal charge>:<id>]<bond>" for monomer atoms and "<bond>[*:<id>]" for adjacent neighboring atoms (NOTE: This functionality is experimental!)

>>> PE_substructs = {
...     "PE": [
...         "[#6D4+0:2](-[#1D1+0:3])(-[#1D1+0:4])(-[#6D4+0:5](-[#1D1+0:6])(-[#1D1+0:7])-[*:8])-[*:1]",
...         "[#6D4+0:2](-[#1D1+0:3])(-[#1D1+0:4])(-[#6D4+0:5](-[#1D1+0:6])(-[#1D1+0:7])-[#1D1+0:8])-[*:1]",
...         "[#6D4+0:2](-[#1D1+0:3])(-[#1D1+0:4])(-[#6D4+0:5](-[#1D1+0:6])(-[#1D1+0:7])-[*:8])-[#1D1+0:1]",
...     ]
... }
>>> top = Topology.from_pdb(
...          get_data_file_path("systems/test_systems/PE.pdb"),
...          _custom_substructures=PE_substructs,
...      )
to_openmm(ensure_unique_atom_names: str | bool = 'residues') openmm.app.Topology[source]

Create an OpenMM Topology object.

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

Contiguously-indexed atoms with the same residue_name, residue_number, insertion_code, 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

Whether to generate new atom names to ensure uniqueness within a molecule or hierarchy element.

  • If the name of a HierarchyScheme is given as a string, new atom names will be generated so that each element of that scheme has unique atom names. Molecules without the given hierarchy scheme will be given unique atom names within that molecule.

  • If True, new atom names will be generated so that atom names are unique within a molecule.

  • If False, the existing atom names will be used.

Returns:

openmm_topology – An OpenMM Topology object

to_file(file: Path | str | TextIO, positions: OMMQuantity | Quantity | ndarray[Any, dtype[_ScalarType_co]] | None = None, file_format: Literal['PDB'] = 'PDB', keep_ids: bool = False, ensure_unique_atom_names: str | bool = 'residues')[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:
  • file – A file-like object to write to, or a path to save the file to.

  • positions

    May 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 dimensions of length

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

    • None (the default), in which case the first conformer of each molecule in the topology will be used.

  • file_format – Output file format. Case insensitive. Currently only supported values are "PDB".

  • keep_ids – If True, keep the residue and chain IDs specified in the Topology rather than generating new ones.

  • ensure_unique_atom_names

    Whether to generate new atom names to ensure uniqueness within a molecule or hierarchy element.

    • If the name of a HierarchyScheme is given as a string, new atom names will be generated so that each element of that scheme has unique atom names. Molecules without the given hierarchy scheme will be given unique atom names within that molecule.

    • If True, new atom names will be generated so that atom names are unique within a molecule.

    • If False, the existing atom names will be used.

    Note that this option cannot guarantee name uniqueness for formats like PDB that truncate long atom names.

get_positions() Quantity | None[source]

Copy the positions of the topology into a new array.

Topology positions are stored as the first conformer of each molecule. If any molecule has no conformers, this method returns None. Note that modifying the returned array will not update the positions in the topology. To change the positions, use Topology.set_positions().

clear_positions()[source]

Clear the positions of this topology by removing all conformers from its consituent molecules.

Note that all conformers will be deleted (in-place) from all molecules. Use Topology.get_positions() if you wish to save them before clearing.

set_positions(array: Quantity)[source]

Set the positions in a topology by copying from a single n×3 array.

Note that modifying the original array will not update the positions in the topology; it must be passed again to set_positions().

Parameters:

array – Positions for the topology. Should be a unit-wrapped array-like object with shape (n_atoms, 3) and dimensions of length.

classmethod from_mdtraj(mdtraj_topology: Topology, unique_molecules: Iterable[Molecule | FrozenMolecule | _SimpleMolecule] | None = None, positions: None | OMMQuantity | Quantity = None)[source]

Construct an OpenFF Topology from an MDTraj Topology

This method guarantees that the order of atoms in the input MDTraj 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.

Hierarchy schemes are taken from the MDTraj topology, not from unique_molecules.

Parameters:
  • mdtraj_topology – The MDTraj Topology object to convert

  • unique_molecules – An iterable containing all the unique molecules in the topology. This is used to identify the molecules in the MDTraj topology and provide any missing chemical information. Each chemical species in the topology must be specified exactly once, though the topology may have any number of copies, including zero. The chemical elements of atoms and their bond connectivity will be used to match these reference molecules to the molecules appearing in the topology. If bond orders are specified in the topology, these will be used in matching as well.

  • positions – Positions for the atoms in the new topology.

Returns:

topology – An OpenFF Topology object, constructed from the molecules in unique_molecules, with the same atom order as the input topology.

Raises:
get_bond_between(i: int | Atom, j: int | Atom) Bond[source]

Returns the bond between two atoms

Parameters:
  • i – Atoms or atom indices to check

  • j – Atoms or atom indices to check

Returns:

bond – The bond between i and j.

is_bonded(i: int | Atom, j: int | Atom) bool[source]

Returns True if the two atoms are bonded

Parameters:
  • i – Atoms or atom indices to check

  • j – Atoms or atom indices to check

Returns:

is_bonded – True if atoms are bonded, False otherwise.

atom(atom_topology_index: int) Atom[source]

Get the Atom at a given Topology atom index.

Parameters:

atom_topology_index – The index of the Atom in this Topology

Returns:

An openff.toolkit.topology.Atom

bond(bond_topology_index: int) Bond[source]

Get the Bond at a given Topology bond index.

Parameters:

bond_topology_index – The index of the Bond in this Topology

Returns:

An openff.toolkit.topology.Bond

add_molecule(molecule: Molecule | FrozenMolecule | _SimpleMolecule) int[source]

Add a copy of the molecule to the topology

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 – Atoms to mark as constrained These atoms may be bonded or not in the Topology

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

  • distance – 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 – Indices of atoms to mark as constrained.

  • jatom – Indices of atoms to mark as constrained.

Returns:

distance – True if constrained but constraints have not yet been applied Distance if constraint has already been added to Topology

classmethod from_bson(serialized)

Instantiate an object from a BSON serialized representation.

Specification: http://bsonspec.org/

Parameters:

serialized – A BSON serialized representation of the object

Returns:

instance – An instantiated object

classmethod from_json(serialized: str)

Instantiate an object from a JSON serialized representation.

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

Parameters:

serialized – A JSON serialized representation of the object

Returns:

instance – An instantiated object

classmethod from_messagepack(serialized)

Instantiate an object from a MessagePack serialized representation.

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

Parameters:

serialized – A MessagePack-encoded bytes serialized representation

Returns:

instance – 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 – A pickled representation of the object

Returns:

instance – An instantiated object

classmethod from_toml(serialized)

Instantiate an object from a TOML serialized representation.

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

Parameters:

serlialized – A TOML serialized representation of the object

Returns:

instance – An instantiated object

classmethod from_xml(serialized)

Instantiate an object from an XML serialized representation.

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

Parameters:

serialized – An XML serialized representation

Returns:

instance – Instantiated object.

classmethod from_yaml(serialized)

Instantiate from a YAML serialized representation.

Specification: http://yaml.org/

Parameters:

serialized – A YAML serialized representation of the object

Returns:

instance – Instantiated object

to_bson()

Return a BSON serialized representation.

Specification: http://bsonspec.org/

Returns:

serialized – A BSON serialized representation of the objecft

to_json(indent=None) str

Return a JSON serialized representation.

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

Parameters:

indent – If not None, will pretty-print with specified number of spaces for indentation

Returns:

serialized – A JSON serialized representation of the object

to_messagepack()

Return a MessagePack representation.

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

Returns:

serialized – 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 – A pickled representation of the object

to_toml()

Return a TOML serialized representation.

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

Returns:

serialized – A TOML serialized representation of the object

to_xml(indent=2)

Return an XML representation.

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

Parameters:

indent – If not None, will pretty-print with specified number of spaces for indentation

Returns:

serialized – A MessagePack-encoded bytes serialized representation.

to_yaml()

Return a YAML serialized representation.

Specification: http://yaml.org/

Returns:

serialized – A YAML serialized representation of the object

visualize(ensure_correct_connectivity: bool = False) NGLWidget[source]

Visualize with NGLView.

Requires all molecules in this topology have positions.

NGLView is a 3D molecular visualization library for use in Jupyter notebooks. Note that for performance reasons, by default the visualized connectivity is inferred from positions and may not reflect the connectivity in the Topology.

Parameters:

ensure_correct_connectivity – If True, the visualization will be guaranteed to reflect the connectivity in the Topology. Note that this will severely degrade performance, especially for topologies with many atoms.

Examples

Visualize a complex PDB file

>>> from openff.toolkit import Topology
>>> from openff.toolkit.utils.utils import get_data_file_path
>>> pdb_filename = get_data_file_path("systems/test_systems/T4_lysozyme_water_ions.pdb")
>>> topology = Topology.from_pdb(pdb_filename)
>>> topology.visualize()  
hierarchy_iterator(iter_name: str) Iterator[HierarchyElement][source]

Iterate over all molecules with the given hierarchy scheme.

Get an iterator over hierarchy elements from all of the molecules in this topology that provide the appropriately named iterator. This iterator will yield hierarchy elements sorted first by the order that molecules are listed in the Topology, and second by the specific sorting of hierarchy elements defined in each molecule. Molecules without the named iterator are not included.

Parameters:

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

Returns:

iterator of HierarchyElement