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)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 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[, ...])Construct an OpenFF
Topology
from an MDTrajTopology
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_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
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.
particle_index
(particle)Deprecated since version 0.11.0.
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.
Return a MessagePack representation.
to_openmm
([ensure_unique_atom_names])Create an OpenMM Topology object.
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
Iterate over improper torsions in the molecule, but only those with trivalent centers, reporting the central atom first in each improper.
iterator over the angles in this Topology.
Get the aromaticity model applied to all molecules in the topology.
Returns an iterator over the atoms in this Topology.
Returns an iterator over the bonds in this Topology
Return the box vectors of the topology, if specified
Returns the constrained atom pairs of the Topology
Returns groups of chemically identical molecules, identified by index and atom map.
iterator over the possible improper torsions in this Topology.
Return whether or not this Topology is intended to be described with periodic boundary conditions.
Returns an iterator over all the Molecules in this Topology
number of angles in this Topology.
Returns the number of atoms in in this Topology.
Returns the number of Bonds in in this Topology.
number of possible improper torsions in this Topology.
Returns the number of molecules in this Topology
Deprecated since version 0.11.0.
number of proper torsions in this Topology.
Deprecated since version 0.11.0.
Deprecated since version 0.11.0.
Deprecated since version 0.11.0.
Deprecated since version 0.11.0.
Deprecated since version 0.11.0.
Returns the number of unique molecules in this Topology
Deprecated since version 0.11.0.
iterator over the proper torsions in this Topology.
Deprecated since version 0.11.0.
Iterate over improper torsions in the molecule, but only those with trivalent centers, reporting the central atom second in each improper.
Deprecated since version 0.11.0.
Deprecated since version 0.11.0.
Deprecated since version 0.11.0.
Deprecated since version 0.11.0.
Get a list of chemically unique molecules in this Topology.
- property unique_molecules: Iterator[Molecule]
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: 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.
- 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[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: int
Returns the number of molecules in this Topology
- 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 angles
iterator over the angles in this Topology.
- Type
Iterable of Tuple[Atom]
- property propers
iterator over the proper torsions in this Topology.
- Type
Iterable of Tuple[TopologyAtom]
- 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.
See also
- 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.
See also
- 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=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
- Returns
matches (list of Topology._ChemicalEnvironmentMatch) – 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 (dict) – A dictionary representation of the topology.
- Returns
topology (Topology) – A Topology created from the dictionary representation
- classmethod from_openmm(openmm_topology: openmm.app.Topology, unique_molecules: Optional[Iterable[FrozenMolecule]] = None, positions: Union[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
.- 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
MissingUniqueMoleculesError – If
unique_molecules
isNone
DuplicateUniqueMoleculeError – If the same connectivity graph is represented by two different molecules in
unique_molecules
ValueError – If a chemically impossible molecule is detected in the topology
- to_openmm(ensure_unique_atom_names: Union[str, bool] = 'residues')[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 (openmm.app.Topology) – An OpenMM Topology object
- to_file(file: Union[Path, str, TextIO], positions: Optional[Union[OMMQuantity, Quantity, ndarray[Any, dtype[ScalarType]]]] = None, file_format: Literal['PDB'] = 'PDB', keep_ids: bool = False, ensure_unique_atom_names: Union[str, bool] = 'residues')[source]
Save coordinates and topology to a PDB file.
Reference: https://github.com/openforcefield/openff-toolkit/issues/502
Notes:
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.
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 (Array with shape
(n_atoms, 3)
and dimensions of length) –May be a…
openmm.unit.Quantity
object which has atomic positions as a list of unit-taggedVec3
objectsopenff.units.unit.Quantity
object which wraps anumpy.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 value is
"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() Optional[Quantity] [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, useTopology.set_positions()
.See also
- 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.
See also
- classmethod from_mdtraj(mdtraj_topology: Topology, unique_molecules: Optional[Iterable[FrozenMolecule]] = None, positions: Union[None, OMMQuantity, Quantity] = None)[source]
Construct an OpenFF
Topology
from an MDTrajTopology
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
MissingUniqueMoleculesError – If
unique_molecules
isNone
DuplicateUniqueMoleculeError – If the same connectivity graph is represented by two different molecules in
unique_molecules
ValueError – If a chemically impossible molecule is detected in the topology
- get_bond_between(i, j)[source]
Returns the bond between two atoms
- is_bonded(i, j)[source]
Returns True if the two atoms are bonded
- 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_molecule(molecule: Union[Molecule, _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 (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 determinedFalse
if constraint is to be removed
- is_constrained(iatom, jatom)[source]
Check if a pair of atoms are marked as constrained.
- 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
- property n_topology_atoms: int
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.n_atoms()
instead.
- property topology_atoms
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.atoms()
instead.
- property n_topology_bonds: int
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.n_bonds()
instead.
- property topology_bonds
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.bonds()
instead.
- property n_topology_particles: int
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.n_particles()
instead.
- property topology_particles
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.particles()
instead.
- property reference_molecules: Iterator[Molecule]
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.unique_molecules()
instead.
- property n_reference_molecules: int
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.n_unique_molecules()
instead.
- property n_topology_molecules: int
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.n_molecules()
instead.
- property topology_molecules
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.molecules()
instead.
- 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
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.n_atoms()
instead.
- 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
Deprecated since version 0.11.0: This property has been deprecated and will soon be removed. Use
Topology.atoms()
instead.
- particle_index(particle) int [source]
Deprecated since version 0.11.0: This method has been deprecated and will soon be removed. Use
Topology.atom_index()
instead.