Interchange

pydantic model openff.interchange.components.interchange.Interchange[source]

Bases: DefaultModel

A object for storing, manipulating, and converting molecular mechanics data.

Warning

This object is in an early and experimental state and unsuitable for production.

Warning

This API is experimental and subject to change.

Fields
Validators
field collections: dict[str, Collection] = {}
field topology: Topology = None
Validated by
field mdconfig: MDConfig = None
field box: openff.models.types.ArrayQuantity = None
Validated by
field positions: openff.models.types.ArrayQuantity = None
field velocities: openff.models.types.ArrayQuantity = None
validator validate_box  »  box[source]
validator validate_topology  »  topology[source]
classmethod from_smirnoff(force_field: ForceField, topology: Topology | list[Molecule], box=None, positions=None, charge_from_molecules: Optional[list[Molecule]] = None, partial_bond_orders_from_molecules: Optional[list[Molecule]] = None, allow_nonintegral_charges: bool = False) Interchange[source]

Create a new object by parameterizing a topology with a SMIRNOFF force field.

Parameters
  • force_field (openff.toolkit.ForceField) – The force field to parameterize the topology with.

  • topology (openff.toolkit.Topology or List[openff.toolkit.Molecule]) – The topology to parameterize, or a list of molecules to construct a topology from and parameterize.

  • box (openff.units.Quantity, optional) – The box vectors associated with the Interchange. If None, box vectors are taken from the topology, if present.

  • positions (openff.units.Quantity, optional) – The positions associated with atoms in the input topology. If None, positions are taken from the molecules in topology, if present on all molecules.

  • charge_from_molecules (List[openff.toolkit.molecule.Molecule], optional) – If specified, partial charges will be taken from the given molecules instead of being determined by the force field.

  • partial_bond_orders_from_molecules (List[openff.toolkit.molecule.Molecule], optional) – If specified, partial bond orders will be taken from the given molecules instead of being determined by the force field.

  • allow_nonintegral_charges (bool, optional, default=False) – If True, allow molecules to have approximately non-integral charges.

Notes

If the Molecule objects in the topology argument each contain conformers, the returned Interchange object will have its positions set via concatenating the 0th conformer of each Molecule.

Examples

Generate an Interchange object from a single-molecule (OpenFF) topology and OpenFF 2.0.0 “Sage”

>>> from openff.interchange import Interchange
>>> from openff.toolkit import ForceField, Molecule
>>> mol = Molecule.from_smiles("CC")
>>> mol.generate_conformers(n_conformers=1)
>>> sage = ForceField("openff-2.0.0.offxml")
>>> interchange = Interchange.from_smirnoff(topology=[mol], force_field=sage)
>>> interchange
Interchange with 8 atoms, non-periodic topology
visualize(backend: str = 'nglview', include_virtual_sites: bool = False) NGLWidget[source]

Visualize this Interchange.

This currently only uses NGLview. Other engines may be added in the future.

Parameters
  • backend (str, default=”nglview”) – The backend to use for visualization. Currently only “nglview” is supported.

  • include_virtual_sites (bool, default=False) – Whether or not to include virtual sites in the visualization.

Returns

widget – The NGLWidget containing the visualization.

Return type

nglview.NGLWidget

minimize(engine: str = 'openmm', force_tolerance: Quantity = _DEFAULT_ENERGY_MINIMIZATION_TOLERANCE, max_iterations: int = 10_000)[source]

Minimize the energy of the system using an available engine.

Updates positions in-place.

Parameters
  • engine (str, default=”openmm”) – The engine to use for minimization. Currently only “openmm” is supported.

  • force_tolerance (openff.units.Quantity, default=10.0 kJ / mol / nm) – The force tolerance to run until during energy minimization.

  • max_iterations (int, default=10_000) – The maximum number of iterations to run during energy minimization.

to_gromacs(prefix: str, decimal: int = 3, hydrogen_mass: float = 1.007947)[source]

Export this Interchange object to GROMACS files.

Parameters
  • prefix (str) – The prefix to use for the GROMACS topology and coordinate files, i.e. “foo” will produce “foo.top” and “foo.gro”.

  • decimal (int, default=3) – The number of decimal places to use when writing the GROMACS coordinate file.

  • hydrogen_mass (float, default=1.007947) – The mass to use for hydrogen atoms if not present in the topology. If non-trivially different than the default value, mass will be transferred from neighboring heavy atoms. Note that this is currently not applied to any waters and is unsupported when virtual sites are present.

to_top(file_path: Path | str, hydrogen_mass: float = 1.007947)[source]

Export this Interchange to a GROMACS topology file.

Parameters
  • file_path – The path to the GROMACS topology file to write.

  • hydrogen_mass (float, default=1.007947) – The mass to use for hydrogen atoms if not present in the topology. If non-trivially different than the default value, mass will be transferred from neighboring heavy atoms. Note that this is currently not applied to any waters and is unsupported when virtual sites are present.

to_gro(file_path: Path | str, decimal: int = 3)[source]

Export this Interchange object to a GROMACS coordinate file.

Parameters
  • file_path (Union[Path, str]) – The path to the GROMACS coordinate file to write.

  • decimal (int, default=3) – The number of decimal places to use when writing the GROMACS coordinate file.

to_lammps(file_path: Path | str, writer='internal')[source]

Export this Interchange to a LAMMPS data file.

to_openmm_system(combine_nonbonded_forces: bool = True, add_constrained_forces: bool = False, ewald_tolerance: float = 1e-4, hydrogen_mass: float = 1.007947)[source]

Export this Interchange to an OpenMM System.

Parameters
  • combine_nonbonded_forces (bool, default=True) – If True, an attempt will be made to combine all non-bonded interactions into a single openmm.NonbondedForce. If False, non-bonded interactions will be split across multiple forces.

  • add_constrained_forces (bool, default=False,) – If True, add valence forces that might be overridden by constraints, i.e. call addBond or addAngle on a bond or angle that is fully constrained.

  • ewald_tolerance (float, default=1e-4) – The value passed to NonbondedForce.setEwaldErrorTolerance

  • hydrogen_mass (float, default=1.007947) – The mass to use for hydrogen atoms if not present in the topology. If non-trivially different than the default value, mass will be transferred from neighboring heavy atoms. Note that this is currently not applied to any waters and is unsupported when virtual sites are present.

Returns

system – The OpenMM System object.

Return type

openmm.System

to_openmm(combine_nonbonded_forces: bool = True, add_constrained_forces: bool = False, ewald_tolerance: float = 1e-4, hydrogen_mass: float = 1.007947)

Export this Interchange to an OpenMM System.

Parameters
  • combine_nonbonded_forces (bool, default=True) – If True, an attempt will be made to combine all non-bonded interactions into a single openmm.NonbondedForce. If False, non-bonded interactions will be split across multiple forces.

  • add_constrained_forces (bool, default=False,) – If True, add valence forces that might be overridden by constraints, i.e. call addBond or addAngle on a bond or angle that is fully constrained.

  • ewald_tolerance (float, default=1e-4) – The value passed to NonbondedForce.setEwaldErrorTolerance

  • hydrogen_mass (float, default=1.007947) – The mass to use for hydrogen atoms if not present in the topology. If non-trivially different than the default value, mass will be transferred from neighboring heavy atoms. Note that this is currently not applied to any waters and is unsupported when virtual sites are present.

Returns

system – The OpenMM System object.

Return type

openmm.System

to_openmm_topology(ensure_unique_atom_names: str | bool = 'residues')[source]

Export components of this Interchange to an OpenMM Topology.

to_openmm_simulation(integrator: openmm.Integrator, combine_nonbonded_forces: bool = True, add_constrained_forces: bool = False, additional_forces: Iterable[openmm.Force] = tuple(), **kwargs) Simulation[source]

Export this Interchange to an OpenMM Simulation object.

Positions are set on the Simulation if present on the Interchange.

Parameters
  • integrator (subclass of openmm.Integrator) – The integrator to use for the simulation.

  • combine_nonbonded_forces (bool, default=False) – If True, an attempt will be made to combine all non-bonded interactions into a single openmm.NonbondedForce. If False, non-bonded interactions will be split across multiple forces.

  • add_constrained_forces (bool, default=False,) – If True, add valence forces that might be overridden by constraints, i.e. call addBond or addAngle on a bond or angle that is fully constrained.

  • additional_forces (Iterable[openmm.Force], default=tuple()) – Additional forces to be added to the system, i.e. barostats that are not added by the force field.

  • **kwargs – Further keyword parameters are passed on to Simulation.__init__()

Returns

simulation – The OpenMM simulation object, possibly with positions set.

Return type

openmm.app.Simulation

Examples

Create an OpenMM simulation with a Langevin integrator:

>>> import openmm
>>> import openmm.unit
>>>
>>> integrator = openmm.LangevinMiddleIntegrator(
...     293.15 * openmm.unit.kelvin,
...     1.0 / openmm.unit.picosecond,
...     2.0 * openmm.unit.femtosecond,
... )
>>> barostat = openmm.MonteCarloBarostat(
...     1.00 * openmm.unit.bar,
...     293.15 * openmm.unit.kelvin,
...     25,
... )
>>> simulation = interchange.to_openmm_simulation(
...     integrator=integrator,
...     additional_forces=[barostat],
... )

Add a barostat:

Re-initializing the Context after adding a Force is necessary due to implementation details in OpenMM. For more, see https://github.com/openmm/openmm/wiki/Frequently-Asked-Questions#why-does-it-ignore-changes-i-make-to-a-system-or-force

to_prmtop(file_path: Path | str, writer='internal')[source]

Export this Interchange to an Amber .prmtop file.

to_pdb(file_path: Path | str, include_virtual_sites: bool = False)[source]

Export this Interchange to a .pdb file.

to_psf(file_path: Path | str)[source]

Export this Interchange to a CHARMM-style .psf file.

to_crd(file_path: Path | str)[source]

Export this Interchange to a CHARMM-style .crd file.

to_inpcrd(file_path: Path | str, writer='internal')[source]

Export this Interchange to an Amber .inpcrd file.

classmethod from_foyer(force_field: Forcefield, topology: Topology, box=None, positions=None, **kwargs) Interchange[source]

Create an Interchange object from a Foyer force field and an OpenFF topology.

Examples

Generate an Interchange object from a single-molecule (OpenFF) topology and the Foyer implementation of OPLS-AA

>>> from openff.interchange import Interchange
>>> from openff.toolkit import Molecule, Topology
>>> from foyer import Forcefield
>>> mol = Molecule.from_smiles("CC")
>>> mol.generate_conformers(n_conformers=1)
>>> top = Topology.from_molecules([mol])
>>> oplsaa = Forcefield(name="oplsaa")
>>> interchange = Interchange.from_foyer(topology=top, force_field=oplsaa)
>>> interchange
Interchange with 8 atoms, non-periodic topology
classmethod from_gromacs(topology_file: Path | str, gro_file: Path | str) Interchange[source]

🧪 Create an Interchange object from GROMACS files.

Experimental

This object is experimental and should not be used in production.

WARNING! This method is experimental and not suitable for production.

Parameters
  • topology_file (Union[Path, str]) – The path to a GROMACS topology file.

  • gro_file (Union[Path, str]) – The path to a GROMACS coordinate file.

Returns

interchange – An Interchange object representing the contents of the GROMACS files.

Return type

Interchange

classmethod from_openmm(system: openmm.System, topology: Optional[Union[openmm.app.Topology, Topology]] = None, positions: pint.Quantity | None = None, box_vectors: pint.Quantity | None = None) Interchange[source]

🧪 Create an Interchange object from OpenMM objects.

Experimental

This object is experimental and should not be used in production.

WARNING! This method is experimental and not yet suitable for production.

Notes

If (topological) bonds in water are missing (physics) parameters, as is often the case with rigid water, these parameters will be filled in with values of 1 Angstrom equilibrium bond length and a default force constant of 50,000 kcal/mol/A^2, representing an arbitrarily stiff harmonic bond, and angle parameters of 104.5 degrees and 1.0 kcal/mol/rad^2, representing an arbitrarily harmonic angle. It is expected that these values will be overwritten by runtime MD options.

Parameters
  • system (openmm.System) – The OpenMM system.

  • topology (openmm.app.Topology, optional) – The OpenMM topology.

  • positions (openmm.unit.Quantity or openff.units.Quantity, optional) – The positions of particles in this system and/or topology.

  • box_vectors (openmm.unit.Quantity or openff.units.Quantity, optional) – The vectors of the simulation box associated with this system and/or topology.

Returns

interchange – An Interchange object representing the contents of the OpenMM objects.

Return type

Interchange

combine(other: Interchange) Interchange[source]

🧪 Combine two Interchange objects. This method is unstable and not yet unsafe.

Experimental

This object is experimental and should not be used in production.