openff.interchange.components.interchange.Interchange
- pydantic model openff.interchange.components.interchange.Interchange[source]
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.
- Config
arbitrary_types_allowed: bool = True
json_encoders: dict = {<class ‘openff.units.units.Quantity’>: <function custom_quantity_encoder at 0x7fdff3466f70>}
json_loads: function = <function json_loader at 0x7fdff3459280>
validate_assignment: bool = True
- classmethod from_smirnoff(force_field: openff.toolkit.typing.engines.smirnoff.forcefield.ForceField, topology: openff.interchange.components.mdtraj._OFFBioTop, box=None) openff.interchange.components.interchange.Interchange [source]
Create a new object by parameterizing a topology with a SMIRNOFF force field.
- Parameters
force_field – The force field to parameterize the topology with.
topology – The topology to parameterize.
box – The box vectors associated with the interchange.
Examples
Generate an Interchange object from a single-molecule (OpenFF) topology and OpenFF 1.0.0 “Parsley”
>>> from openff.interchange.components.interchange import Interchange >>> from openff.interchange.components.mdtraj import _OFFBioTop >>> from openff.toolkit.topology import Molecule >>> from openff.toolkit.typing.engines.smirnoff import ForceField >>> import mdtraj as md >>> mol = Molecule.from_smiles("CC") >>> mol.generate_conformers(n_conformers=1) >>> top = _OFFBioTop.from_molecules([mol]) >>> top.mdtop = md.Topology.from_openmm(top.to_openmm()) >>> parsley = ForceField("openff-1.0.0.offxml") >>> interchange = Interchange.from_smirnoff(topology=top, force_field=parsley) >>> interchange Interchange with 8 atoms, non-periodic topology
- to_gro(file_path: Union[pathlib.Path, str], writer='internal', decimal: int = 8)[source]
Export this Interchange object to a .gro file.
- to_top(file_path: Union[pathlib.Path, str], writer='internal')[source]
Export this interchange to a .top file.
- to_lammps(file_path: Union[pathlib.Path, str], writer='internal')[source]
Export this Interchange to a LAMMPS data file.
- to_openmm(combine_nonbonded_forces: bool = False)[source]
Export this interchange to an OpenMM System.
- to_prmtop(file_path: Union[pathlib.Path, str], writer='internal')[source]
Export this interchange to an Amber .prmtop file.
- to_psf(file_path: Union[pathlib.Path, str])[source]
Export this interchange to a CHARMM-style .psf file.
- to_crd(file_path: Union[pathlib.Path, str])[source]
Export this interchange to a CHARMM-style .crd file.
- to_inpcrd(file_path: Union[pathlib.Path, str], writer='internal')[source]
Export this interchange to an Amber .inpcrd file.
- classmethod from_foyer(topology: _OFFBioTop, force_field: FoyerForcefield, **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.components.interchange import Interchange >>> from openff.interchange.components.mdtraj import _OFFBioTop >>> from openff.toolkit.topology import Molecule >>> from foyer import Forcefield >>> import mdtraj as md >>> mol = Molecule.from_smiles("CC") >>> mol.generate_conformers(n_conformers=1) >>> top = _OFFBioTop.from_molecules([mol]) >>> top.mdtop = md.Topology.from_openmm(top.to_openmm()) >>> 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: Union[pathlib.Path, str], gro_file: Union[pathlib.Path, str], reader='intermol') openff.interchange.components.interchange.Interchange [source]
Create an Interchange object from GROMACS files.