openforcefield.typing.engines.smirnoff.forcefield.
ForceField
(*sources, parameter_handler_classes=None, parameter_io_handler_classes=None, disable_version_check=False)[source]¶A factory that assigns SMIRNOFF parameters to a molecular system
ForceField
is a factory that constructs an OpenMM simtk.openmm.System
object from a
openforcefield.topology.Topology
object defining a (bio)molecular system containing one or more molecules.
When a ForceField
object is created from one or more specified SMIRNOFF serialized representations,
all ParameterHandler
subclasses currently imported are identified and registered to handle different
sections of the SMIRNOFF force field definition file(s).
All ParameterIOHandler
subclasses currently imported are identified and registered to handle different
serialization formats (such as XML).
The force field definition is processed by these handlers to populate the ForceField
object model data
structures that can easily be manipulated via the API:
Processing a Topology
object defining a chemical system will then call all :class`ParameterHandler`
objects in an order guaranteed to satisfy the declared processing order constraints of each
:class`ParameterHandler`.
Examples
Create a new ForceField containing the smirnoff99Frosst parameter set:
>>> from openforcefield.typing.engines.smirnoff import ForceField
>>> forcefield = ForceField('smirnoff99Frosst.offxml')
Create an OpenMM system from a openforcefield.topology.Topology
object:
>>> from openforcefield.topology import Molecule, Topology
>>> ethanol = Molecule.from_smiles('CCO')
>>> topology = Topology.from_molecules(molecules=[ethanol])
>>> system = forcefield.create_openmm_system(topology)
Modify the long-range electrostatics method:
>>> forcefield.get_handler('Electrostatics').method = 'PME'
Inspect the first few vdW parameters:
>>> low_precedence_parameters = forcefield.get_handler('vdW').parameters[0:3]
Retrieve the vdW parameters by SMIRKS string and manipulate it:
>>> parameter = forcefield.get_handler('vdW').parameters['[#1:1]-[#7]']
>>> parameter.sigma += 0.1 * unit.angstroms
>>> parameter.epsilon *= 1.02
Make a child vdW type more specific (checking modified SMIRKS for validity):
>>> forcefield.get_handler('vdW').parameters[-1].smirks += '~[#53]'
Warning
While we check whether the modified SMIRKS is still valid and has the appropriate valence type, we currently don’t check whether the typing remains hierarchical, which could result in some types no longer being assignable because more general types now come below them and preferentially match.
Delete a parameter:
>>> del forcefield.get_handler('vdW').parameters['[#1:1]-[#6X4]']
Insert a parameter at a specific point in the parameter tree:
>>> from openforcefield.typing.engines.smirnoff import vdWHandler
>>> new_parameter = vdWHandler.vdWType(smirks='[*:1]', epsilon=0.0157*unit.kilocalories_per_mole, rmin_half=0.6000*unit.angstroms)
>>> forcefield.get_handler('vdW').parameters.insert(0, new_parameter)
Warning
We currently don’t check whether removing a parameter could accidentally remove the root type, so it’s possible to no longer type all molecules this way.
Attributes: |
|
---|
Methods
create_openmm_system (topology, **kwargs) |
Create an OpenMM System representing the interactions for the specified Topology with the current force field |
create_parmed_structure (topology, positions, …) |
Create a ParmEd Structure object representing the interactions for the specified Topology with the current force field |
get_handler (tagname[, handler_kwargs]) |
Retrieve the parameter handlers associated with the provided tagname. |
get_io_handler (io_format) |
Retrieve the parameter handlers associated with the provided tagname. |
label_molecules (topology) |
Return labels for a list of molecules corresponding to parameters from this force field. |
load_smirnoff_data (smirnoff_data) |
Add parameters from a SMIRNOFF-format data structure to this ForceField. |
parse_smirnoff_from_source (source) |
Reads a SMIRNOFF data structure from a source, which can be one of many types. |
parse_sources (sources) |
Parse a SMIRNOFF force field definition. |
register_parameter_handler (…) |
Register a new ParameterHandler from a specified class, instantiating the ParameterHandler object and making it available for lookup in the ForceField. |
register_parameter_io_handler (…) |
Register a new ParameterIOHandler from a specified class, instantiating the ParameterIOHandler object and making it available for lookup in the ForceField. |
to_smirnoff_data () |
Convert this ForceField and all related ParameterHandlers to an OrderedDict representing a SMIRNOFF data object. |
__init__
(*sources, parameter_handler_classes=None, parameter_io_handler_classes=None, disable_version_check=False)[source]¶Create a new ForceField
object from one or more SMIRNOFF parameter definition files.
Parameters: |
|
---|
Examples
Load one SMIRNOFF parameter set in XML format (searching the package data directory by default, which includes some standard parameter sets):
>>> forcefield = ForceField('smirnoff99Frosst.offxml')
Load multiple SMIRNOFF parameter sets:
forcefield = ForceField(‘smirnoff99Frosst.offxml’, ‘tip3p.offxml’)
Load a parameter set from a string:
>>> offxml = '<SMIRNOFF version="1.0" aromaticity_model="OEAroModel_MDL"/>'
>>> forcefield = ForceField(offxml)
Methods
__init__ (*sources[, …]) |
Create a new ForceField object from one or more SMIRNOFF parameter definition files. |
create_openmm_system (topology, **kwargs) |
Create an OpenMM System representing the interactions for the specified Topology with the current force field |
create_parmed_structure (topology, positions, …) |
Create a ParmEd Structure object representing the interactions for the specified Topology with the current force field |
get_handler (tagname[, handler_kwargs]) |
Retrieve the parameter handlers associated with the provided tagname. |
get_io_handler (io_format) |
Retrieve the parameter handlers associated with the provided tagname. |
label_molecules (topology) |
Return labels for a list of molecules corresponding to parameters from this force field. |
load_smirnoff_data (smirnoff_data) |
Add parameters from a SMIRNOFF-format data structure to this ForceField. |
parse_smirnoff_from_source (source) |
Reads a SMIRNOFF data structure from a source, which can be one of many types. |
parse_sources (sources) |
Parse a SMIRNOFF force field definition. |
register_parameter_handler (…) |
Register a new ParameterHandler from a specified class, instantiating the ParameterHandler object and making it available for lookup in the ForceField. |
register_parameter_io_handler (…) |
Register a new ParameterIOHandler from a specified class, instantiating the ParameterIOHandler object and making it available for lookup in the ForceField. |
to_smirnoff_data () |
Convert this ForceField and all related ParameterHandlers to an OrderedDict representing a SMIRNOFF data object. |
register_parameter_handler
(parameter_handler_class, parameter_handler_kwargs)[source]¶Register a new ParameterHandler from a specified class, instantiating the ParameterHandler object and making it available for lookup in the ForceField.
Warning
This API is experimental and subject to change.
Parameters: |
|
---|---|
Returns: |
|
register_parameter_io_handler
(parameter_io_handler_class)[source]¶Register a new ParameterIOHandler from a specified class, instantiating the ParameterIOHandler object and making it available for lookup in the ForceField.
Warning
This API is experimental and subject to change.
Parameters: |
|
---|
get_handler
(tagname, handler_kwargs=None)[source]¶Retrieve the parameter handlers associated with the provided tagname.
If the parameter handler has not yet been instantiated, it will be created. If a parameter handler object already exists, it will be checked for compatibility and an Exception raised if it is incompatible with the provided kwargs.
Parameters: |
|
---|---|
Returns: |
|
Raises: |
|
get_io_handler
(io_format)[source]¶Retrieve the parameter handlers associated with the provided tagname.
If the parameter handler has not yet been instantiated, it will be created. If a parameter handler object already exists, it will be checked for compatibility and an Exception raised if it is incompatible with the provided kwargs.
Parameters: |
|
---|---|
Returns: |
|
Raises: |
|
parse_sources
(sources)[source]¶Parse a SMIRNOFF force field definition.
Parameters: |
|
---|
to_smirnoff_data
()[source]¶Convert this ForceField and all related ParameterHandlers to an OrderedDict representing a SMIRNOFF data object.
Returns: |
|
---|
load_smirnoff_data
(smirnoff_data)[source]¶Add parameters from a SMIRNOFF-format data structure to this ForceField.
Parameters: |
|
---|
parse_smirnoff_from_source
(source)[source]¶Reads a SMIRNOFF data structure from a source, which can be one of many types.
Parameters: |
|
---|---|
Returns: |
|
create_openmm_system
(topology, **kwargs)[source]¶Create an OpenMM System representing the interactions for the specified Topology with the current force field
Parameters: |
|
---|---|
Returns: |
|
create_parmed_structure
(topology, positions, **kwargs)[source]¶Create a ParmEd Structure object representing the interactions for the specified Topology with the current force field
This method creates a ParmEd Structure
object containing a topology, positions, and parameters.
Parameters: |
|
---|---|
Returns: |
|
label_molecules
(topology)[source]¶Return labels for a list of molecules corresponding to parameters from this force field. For each molecule, a dictionary of force types is returned, and for each force type, each force term is provided with the atoms involved, the parameter id assigned, and the corresponding SMIRKS.
Parameters: |
|
---|---|
Returns: |
|