ForceField

class openmm.app.forcefield.ForceField(*files)

A ForceField constructs OpenMM System objects based on a Topology.

__init__(*files)

Load one or more XML files and create a ForceField object based on them.

Parameters

files (list) – A list of XML files defining the force field. Each entry may be an absolute file path, a path relative to the current working directory, a path relative to this module’s data subdirectory (for built in force fields), or an open file-like object with a read() method from which the forcefield XML data can be loaded.

Methods

__init__(*files)

Load one or more XML files and create a ForceField object based on them.

createSystem(topology[, nonbondedMethod, …])

Construct an OpenMM System representing a Topology with this force field.

generateTemplatesForUnmatchedResidues(topology)

Generate forcefield residue templates for residues in specified topology for which no forcefield templates are available.

getGenerators()

Get the list of all registered generators.

getMatchingTemplates(topology[, …])

Return a list of forcefield residue templates matching residues in the specified topology.

getUnmatchedResidues(topology[, …])

Return a list of Residue objects from specified topology for which no forcefield templates are available.

loadFile(files[, resname_prefix])

Load an XML file and add the definitions from it to this ForceField.

registerAtomType(parameters)

Register a new atom type.

registerGenerator(generator)

Register a new generator.

registerPatch(patch)

Register a new patch that can be applied to templates.

registerResidueTemplate(template)

Register a new residue template.

registerScript(script)

Register a new script to be executed after building the System.

registerTemplateGenerator(generator)

Register a residue template generator that can be used to parameterize residues that do not match existing forcefield templates.

registerTemplateMatcher(matcher)

Register an object that can override the default logic for matching templates to residues.

registerTemplatePatch(residue, patch, …)

Register that a particular patch can be used with a particular residue.

loadFile(files, resname_prefix='')

Load an XML file and add the definitions from it to this ForceField.

Parameters
  • files (string or file or tuple) – An XML file or tuple of XML files containing force field definitions. Each entry may be either an absolute file path, a path relative to the current working directory, a path relative to this module’s data subdirectory (for built in force fields), or an open file-like object with a read() method from which the forcefield XML data can be loaded.

  • prefix (string) – An optional string to be prepended to each residue name found in the loaded files.

getGenerators()

Get the list of all registered generators.

registerGenerator(generator)

Register a new generator.

registerAtomType(parameters)

Register a new atom type.

registerResidueTemplate(template)

Register a new residue template.

registerPatch(patch)

Register a new patch that can be applied to templates.

registerTemplatePatch(residue, patch, patchResidueIndex)

Register that a particular patch can be used with a particular residue.

registerScript(script)

Register a new script to be executed after building the System.

registerTemplateMatcher(matcher)

Register an object that can override the default logic for matching templates to residues.

A template matcher is a callable object that can be invoked as:

template = f(forcefield, residue, bondedToAtom, ignoreExternalBonds, ignoreExtraParticles)

where forcefield is the ForceField invoking it, residue is an openmm.app.Residue object, bondedToAtom[i] is the set of atoms bonded to atom index i, and ignoreExternalBonds and ignoreExtraParticles indicate whether external bonds and extra particules should be considered in matching templates.

It should return a _TemplateData object that matches the residue. Alternatively it may return None, in which case the standard logic will be used to find a template for the residue.

Caution

This method is experimental, and its API is subject to change.

registerTemplateGenerator(generator)

Register a residue template generator that can be used to parameterize residues that do not match existing forcefield templates.

This functionality can be used to add handlers to parameterize small molecules or unnatural/modified residues.

Caution

This method is experimental, and its API is subject to change.

Parameters
  • generator (function) – A function that will be called when a residue is encountered that does not match an existing forcefield template.

  • a residue without a template is encountered (When) –

  • generator function is called with (the) –

  • :: – success = generator(forcefield, residue)

  • forcefield is the calling ForceField object and residue is a openmm.app.topology.Residue object. (where) –

  • must conform to the following API (generator) –

  • ::

    generator API

    param forcefield

    The ForceField object to which residue templates and/or parameters are to be added.

    type forcefield

    openmm.app.ForceField

    param residue

    The residue topology for which a template is to be generated.

    type residue

    openmm.app.Topology.Residue

    returns
    • success (bool) – If the generator is able to successfully parameterize the residue, True is returned. If the generator cannot parameterize the residue, it should return False and not modify forcefield.

    • The generator should either register a residue template directly with forcefield.registerResidueTemplate(template)

    • or it should call forcefield.loadFile(file) to load residue definitions from an ffxml file.

    • It can also use the ForceField programmatic API to add additional atom types (via forcefield.registerAtomType(parameters))

    • or additional parameters.

getUnmatchedResidues(topology, residueTemplates={})

Return a list of Residue objects from specified topology for which no forcefield templates are available.

Caution

This method is experimental, and its API is subject to change.

Parameters
  • topology (Topology) – The Topology whose residues are to be checked against the forcefield residue templates.

  • residueTemplates (dict=dict()) – Specifies which template to use for particular residues. The keys should be Residue objects from the Topology, and the values should be the names of the templates to use for them. This is useful when a ForceField contains multiple templates that can match the same residue (e.g Fe2+ and Fe3+ templates in the ForceField for a monoatomic iron ion in the Topology).

Returns

  • unmatched_residues (list of Residue) – List of Residue objects from topology for which no forcefield residue templates are available. Note that multiple instances of the same residue appearing at different points in the topology may be returned.

  • This method may be of use in generating missing residue templates or diagnosing parameterization failures.

getMatchingTemplates(topology, ignoreExternalBonds=False)

Return a list of forcefield residue templates matching residues in the specified topology.

Caution

This method is experimental, and its API is subject to change.

Parameters
  • topology (Topology) – The Topology whose residues are to be checked against the forcefield residue templates.

  • ignoreExternalBonds (bool=False) – If true, ignore external bonds when matching residues to templates.

Returns

  • templates (list of _TemplateData) – List of forcefield residue templates corresponding to residues in the topology. templates[index] is template corresponding to residue index in topology.residues()

  • This method may be of use in debugging issues related to parameter assignment.

generateTemplatesForUnmatchedResidues(topology)

Generate forcefield residue templates for residues in specified topology for which no forcefield templates are available.

Caution

This method is experimental, and its API is subject to change.

Parameters

topology (Topology) – The Topology whose residues are to be checked against the forcefield residue templates.

Returns

  • templates (list of _TemplateData) – List of forcefield residue templates corresponding to residues in topology for which no forcefield templates are currently available. Atom types will be set to None, but template name, atom names, elements, and connectivity will be taken from corresponding Residue objects.

  • residues (list of Residue) – List of Residue objects that were used to generate the templates. residues[index] is the Residue that was used to generate the template templates[index]

createSystem(topology, nonbondedMethod=NoCutoff, nonbondedCutoff=Quantity(value=1.0, unit=nanometer), constraints=None, rigidWater=None, removeCMMotion=True, hydrogenMass=None, residueTemplates={}, ignoreExternalBonds=False, switchDistance=None, flexibleConstraints=False, drudeMass=Quantity(value=0.4, unit=dalton), **args)

Construct an OpenMM System representing a Topology with this force field.

Parameters
  • topology (Topology) – The Topology for which to create a System

  • nonbondedMethod (object=NoCutoff) – The method to use for nonbonded interactions. Allowed values are NoCutoff, CutoffNonPeriodic, CutoffPeriodic, Ewald, PME, or LJPME.

  • nonbondedCutoff (distance=1*nanometer) – The cutoff distance to use for nonbonded interactions

  • constraints (object=None) – Specifies which bonds and angles should be implemented with constraints. Allowed values are None, HBonds, AllBonds, or HAngles.

  • rigidWater (boolean=None) – If true, water molecules will be fully rigid regardless of the value passed for the constraints argument. If None (the default), it uses the default behavior for this force field’s water model.

  • removeCMMotion (boolean=True) – If true, a CMMotionRemover will be added to the System

  • hydrogenMass (mass=None) – The mass to use for hydrogen atoms bound to heavy atoms. Any mass added to a hydrogen is subtracted from the heavy atom to keep their total mass the same. If rigidWater is used to make water molecules rigid, then water hydrogens are not altered.

  • residueTemplates (dict=dict()) – Specifies which template to use for particular residues. The keys should be Residue objects from the Topology, and the values should be the names of the templates to use for them. This is useful when a ForceField contains multiple templates that can match the same residue (e.g Fe2+ and Fe3+ templates in the ForceField for a monoatomic iron ion in the Topology).

  • ignoreExternalBonds (boolean=False) – If true, ignore external bonds when matching residues to templates. This is useful when the Topology represents one piece of a larger molecule, so chains are not terminated properly. This option can create ambiguities where multiple templates match the same residue. If that happens, use the residueTemplates argument to specify which one to use.

  • switchDistance (float=None) – The distance at which the potential energy switching function is turned on for Lennard-Jones interactions. If this is None, no switching function will be used.

  • flexibleConstraints (boolean=False) – If True, parameters for constrained degrees of freedom will be added to the System

  • drudeMass (mass=0.4*amu) – The mass to use for Drude particles. Any mass added to a Drude particle is subtracted from its parent atom to keep their total mass the same.

  • args – Arbitrary additional keyword arguments may also be specified. This allows extra parameters to be specified that are specific to particular force fields.

Returns

the newly created System

Return type

system