openforcefield.topology.Topology

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

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

As of the 0.7.0 release, the Topology particle indexing system puts all atoms before all virtualsites. This ensures that atoms keep the same Topology particle index value, even if the Topology is modified during system creation by the addition of virtual sites.

Warning

This API is experimental and subject to change.

Examples

Import some utilities

>>> from simtk.openmm import app
>>> from openforcefield.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 openforcefield.topology 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)]
Attributes
angles

Iterable of Tuple[TopologyAtom]: iterator over the angles in this Topology.

aromaticity_model

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

box_vectors

Return the box vectors of the topology, if specified

charge_model

Get the partial charge model applied to all molecules in the topology.

constrained_atom_pairs

Returns the constrained atom pairs of the Topology

fractional_bond_order_model

Get the fractional bond order model for the Topology.

impropers

Iterable of Tuple[TopologyAtom]: iterator over the improper torsions in this Topology.

n_angles

int: number of angles in this Topology.

n_impropers

int: number of improper torsions in this Topology.

n_propers

int: number of proper torsions in this Topology.

n_reference_molecules

Returns the number of reference (unique) molecules in in this Topology.

n_topology_atoms

Returns the number of topology atoms in in this Topology.

n_topology_bonds

Returns the number of TopologyBonds in in this Topology.

n_topology_molecules

Returns the number of topology molecules in in this Topology.

n_topology_particles

Returns the number of topology particles (TopologyAtoms and TopologyVirtualSites) in in this Topology.

n_topology_virtual_sites

Returns the number of TopologyVirtualSites in in this Topology.

propers

Iterable of Tuple[TopologyAtom]: iterator over the proper torsions in this Topology.

reference_molecules

Get an iterator of reference molecules in this Topology.

topology_atoms

Returns an iterator over the atoms in this Topology.

topology_bonds

Returns an iterator over the bonds in this Topology

topology_molecules

Returns an iterator over all the TopologyMolecules in this Topology

topology_particles

Returns an iterator over the particles (TopologyAtoms and TopologyVirtualSites) in this Topology.

topology_virtual_sites

Get an iterator over the virtual sites in this Topology

Methods

add_constraint(iatom, jatom[, distance])

Mark a pair of atoms as constrained.

add_molecule(molecule[, …])

Add a Molecule to the Topology.

add_particle(particle)

Add a Particle 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 TopologyAtom at a given Topology atom index.

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.

from_bson(serialized)

Instantiate an object from a BSON serialized representation.

from_dict(d)

Static constructor from dictionary representation.

from_json(serialized)

Instantiate an object from a JSON serialized representation.

from_mdtraj(mdtraj_topology[, unique_molecules])

Construct an openforcefield 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 openforcefield Topology object from an OpenMM Topology object.

from_parmed(parmed_structure[, unique_molecules])

Warning

This functionality will be implemented in a future toolkit release.

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

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.

to_bson()

Return a BSON serialized representation.

to_dict()

Convert to dictionary representation.

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_parmed()

Warning

This functionality will be implemented in a future toolkit release.

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.

virtual_site(vsite_topology_index)

Get the TopologyAtom at a given Topology atom index.

__init__(other=None)[source]

Create a new Topology.

Parameters
otheroptional, 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[, …])

Add a Molecule to the Topology.

add_particle(particle)

Add a Particle 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 TopologyAtom at a given Topology atom index.

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.

from_bson(serialized)

Instantiate an object from a BSON serialized representation.

from_dict(d)

Static constructor from dictionary representation.

from_json(serialized)

Instantiate an object from a JSON serialized representation.

from_mdtraj(mdtraj_topology[, unique_molecules])

Construct an openforcefield 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 openforcefield Topology object from an OpenMM Topology object.

from_parmed(parmed_structure[, unique_molecules])

Warning

This functionality will be implemented in a future toolkit release.

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

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.

to_bson()

Return a BSON serialized representation.

to_dict()

Convert to dictionary representation.

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_parmed()

Warning

This functionality will be implemented in a future toolkit release.

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.

virtual_site(vsite_topology_index)

Get the TopologyAtom at a given Topology atom index.

Attributes

angles

Iterable of Tuple[TopologyAtom]: iterator over the angles in this Topology.

aromaticity_model

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

box_vectors

Return the box vectors of the topology, if specified Returns ——- box_vectors : simtk.unit.Quantity wrapped numpy array The unit-wrapped box vectors of this topology

charge_model

Get the partial charge model applied to all molecules in the topology.

constrained_atom_pairs

Returns the constrained atom pairs of the Topology

fractional_bond_order_model

Get the fractional bond order model for the Topology.

impropers

Iterable of Tuple[TopologyAtom]: iterator over the improper torsions in this Topology.

n_angles

int: number of angles in this Topology.

n_impropers

int: number of improper torsions in this Topology.

n_propers

int: number of proper torsions in this Topology.

n_reference_molecules

Returns the number of reference (unique) molecules in in this Topology.

n_topology_atoms

Returns the number of topology atoms in in this Topology.

n_topology_bonds

Returns the number of TopologyBonds in in this Topology.

n_topology_molecules

Returns the number of topology molecules in in this Topology.

n_topology_particles

Returns the number of topology particles (TopologyAtoms and TopologyVirtualSites) in in this Topology.

n_topology_virtual_sites

Returns the number of TopologyVirtualSites in in this Topology.

propers

Iterable of Tuple[TopologyAtom]: iterator over the proper torsions in this Topology.

reference_molecules

Get an iterator of reference molecules in this Topology.

topology_atoms

Returns an iterator over the atoms in this Topology.

topology_bonds

Returns an iterator over the bonds in this Topology

topology_molecules

Returns an iterator over all the TopologyMolecules in this Topology

topology_particles

Returns an iterator over the particles (TopologyAtoms and TopologyVirtualSites) in this Topology.

topology_virtual_sites

Get an iterator over the virtual sites in this Topology