Force fields
The mmelemental.models.forcefield submodule provides models that describe force field definitions or parametrized molecules (e.g. data stored in “topology” files).
ForceField
Description
This is the most basic model for storing force field data that describe atomistic or coarse-grained inter-particle potentials. Model creation occurs with a kwargs constructor as shown by equivalent operations below:
>>> ff = mmelemental.models.ForceField(
**{
"name": "water_ff",
"symbols": ["O", "H", "H"],
"charges": [-0.834, 0.417, 0.417],
"masses": [16.0, 1.008, 1.008],
"defs": ["OW", "HW", "HW"],
}
)
>>> ff
ForceField(name='water_ff', form=[], hash='8352e99')
Note that this force field object has no form since we did not define any bonded (or non-bonded) terms as part of the interaction. This can be specified using additional models explained in the next subsections.
In addition, ForceField provides from_data()
and from_file()
methods to create a forcefield from data (e.g. parmed.Structure) or file objects.
The methods to_data()
and to_file()
enable converting a force field to data and file objects, respectively. See the API for more details.
API
- class mmelemental.models.forcefield.forcefield.ForceField(**kwargs)[source]
- Parameters
schema_name (ConstrainedStrValue, Default: mmschema_forcefield) – The MMSchema specification to which this model conforms. Explicitly fixed as mmschema_forcefield.
schema_version (int, Default: 1) – The version number of
schema_name
to which this model conforms.author (str, Optional) – Author name to assign to the force field this model stores. This field can be arbitrary.
name (str, Optional) – Common or human-readable name to assign to this model. This field can be arbitrary.
version (str, Optional) – Version of the force field this model stores. This field can be arbitrary.
comment (str, Optional) – Additional comments for this model. Intended for pure human/user consumption and clarity.
symbols (Array, Optional) – An ordered (natom,) list of particle (e.g. atomic elemental) symbols.
nonbonded (Union[
NonBonded
,NonBonded
], Optional) – Non-bonded parameters model.bonds (Union[
Bonds
,Bonds
], Optional) – 2-body covalent bond model.angles (Union[
Angles
,Angles
], Optional) – 3-body angular bond model.dihedrals (Union[
Dihedrals
,Dihedrals
], Optional) – 4-body torsional bond model.dihedrals_improper (Union[
DihedralsImproper
,DihedralsImproper
], Optional) – Improper 4-body torsional bond model.charges (Array, Optional) – Atomic charges. Default unit is in elementary charge units.
charges_units (str, Default: e) – Atomic charge unit.
masses (Array, Optional) – List of atomic masses. If not provided, the mass of each atom is inferred from its most common isotope. If this is provided, it must be the same length as
symbols
.masses_units (str, Default: amu) – Units for atomic masses. Defaults to unified atomic mass unit.
charge_groups (Array, Optional) – Charge groups per atom. Length of the array must be natoms.
exclusions (str, Optional) – Which pairs of bonded atoms to exclude from non-bonded calculations. The rules to apply in choosing bonded exclusions are specifed in the configuration file using the exclude parameter. The choices for exclusions are None, 1-2, 1-3, 1-4, etc. With None, no atom pairs are excluded. With 1-2, only atoms that are connected via a linear bond are excluded. With 1-3, any pair of atoms connected via a bond or bonded to a common third atom are excluded. With 1-4, any atoms that are connected by a linear bond, or a sequence of two bonds, or a sequence of three bonds, are excluded. With scaled1-4, exclusions are applied in the same manner as the 1-3 setting, but those pairs that are connected by a sequence of 3 bonds are calculated using the modified 1-4 methods described rather than the standard force calculations.
inclusions (str, Optional) – Which pairs of 1-4 excluded bonded atoms to include in non-bonded calculations.
defs (List[str], Optional) – Particle definition. For atomic forcefields, this could be the atom type (e.g. HH31) or SMIRKS (OFF) representation. The type names are associated with the atomic elements defined in other objects e.g. see the :class:
Molecule
model.substructs (List[Tuple[str, int]], Optional) – A list of substructure names the particles belong to. E.g. [(‘ALA’, 4), (‘ACE’, 0)] means atom1 belong to residue ALA (alanine) with residue number 4, while atom2 belongs to residue ACE (acetyl) with residue number 0.
templates (Dict[List[str]], Optional) – A list of template definitions typically in terms of atom types. E.g. {‘ACE’: [‘HH31’, ‘CH3’, ‘HH32’, ‘HH33’, ‘C’, ‘O’]}.
atomic_numbers (Array, Optional) – An optional ordered 1-D array-like object of atomic numbers of shape (nat,). Index matches the 0-indexed indices of all other per-atom settings like
symbols
. Values are inferred from thesymbols
list if not explicitly set.provenance (
Provenance
, Default: {‘creator’: ‘MMElemental’, ‘version’: ‘v0.1.1’, ‘routine’: ‘mmelemental.models.forcefield.forcefield’}) – The provenance information about how this object (and its attributes) were generated, provided, and manipulated.extras (Dict[Any], Optional) – Additional information to bundle with the object. Use for schema development and scratch space.
- Return type
None
- classmethod from_data(data, **kwargs)[source]
Constructs a ForceField object from a data object.
- Parameters
data (Any) – Data to construct ForceField from.
**kwargs (Optional[Dict[str, Any]], optional) – Additional kwargs to pass to the constructors.
- Returns
A constructed ForceField object.
- Return type
- classmethod from_file(filename, dtype=None, translator=None, **kwargs)[source]
Constructs a ForceField object from a file.
- Parameters
filename (str) – The topology or FF filename to build from.
dtype (Optional[str], optional) – The type of file to interpret e.g. psf. If unset, mmelemental attempts to discover the file type.
translator (Optional[str], optional) – Translator name e.g. mmic_parmed. Takes precedence over dtype. If unset, MMElemental attempts to find an appropriate translator if it is registered in the
TransComponent
class.**kwargs (Optional[Dict[str, Any]], optional) – Any additional keywords to pass to the constructor.
- Returns
A constructed ForceField object.
- Return type
- property is_topology
Returns True if model contains “topological” (i.e. subset of assigned ff params) data rather than forcefield definition.
- to_data(dtype=None, translator=None, **kwargs)[source]
Constructs a toolkit-specific forcefield from MMSchema ForceField. Which toolkit-specific component is called depends on which package is installed on the system.
- Parameters
translator (Optional[str], optional) – Translator name e.g. mmic_parmed. Takes precedence over dtype. If unset, MMElemental attempts to find an appropriate translator if it is registered in the
TransComponent
class.dtype (str, optional) – Data type e.g. mdanalysis, parmed, etc.
**kwargs (Optional[Dict[str, Any]]) – Additional kwargs to pass to the constructors.
Results –
------- –
ToolkitModel – Toolkit-specific ForceField object
kwargs (Dict[str, Any]) –
- Return type
ToolkitModel
- to_file(filename, dtype=None, translator=None, **kwargs)[source]
Writes the ForceField to a file.
- Parameters
filename (str) – The filename to write to
dtype (Optional[str], optional) – The type of file to write (e.g. psf, top, etc.), attempts to infer dtype from file extension if not provided.
translator (Optional[str], optional) – Translator name e.g. mmic_parmed. Takes precedence over dtype. If unset, MMElemental attempts to find an appropriate translator if it is registered in the
TransComponent
class.**kwargs (Optional[str, Dict], optional) – Additional kwargs to pass to the constructor.
kwargs (Dict[str, Any]) –
- Return type
None
NonBonded
Description
The NonBonded model describes potentials for non-connected particles. Model creation occurs with a kwargs constructor as shown by equivalent operations below:
>>> nb = mmelemental.models.forcefield.NonBonded(
**{
"form": "LennardJones",
"params": {"epsilon": [0.636386, 0.0, 0.0], "sigma": [1.575305, 0.0, 0.0]},
}
)
This model can be used as an input to ForceField for an extended description of non-bonded interactions.
>>> ff = mmelemental.models.ForceField(
**{
"name": "water_ff",
"symbols": ["O", "H", "H"],
"charges": [-0.834, 0.417, 0.417],
"masses": [16.0, 1.008, 1.008],
"defs": ["OW", "HW", "HW"],
"nonbonded": nb,
}
)
>>> ff
ForceField(name='water_ff', form=['NonBonded'], hash='0cbf0de')
API
- class mmelemental.models.forcefield.nonbonded.NonBonded(**data)[source]
Model that describes non-bonded interactions between particles.
- Parameters
form (str) – Name or form of the potential. See cls.supported_potentials for available potentials.
version (str, Optional) – Version of the force field this model stores. This field can be arbitrary.
params (Any) – Specific force field parameters model.
defs (List[str], Optional) – Particle definition. For atomic forcefields, this could be the atom type (e.g. HH31) or SMIRKS (OFF) representation.
extras (Dict[Any], Optional) – Additional information to bundle with the object. Use for schema development and scratch space.
combination_rule (str, Default: Lorentz-Berthelot) – Combination rule for the force field.
- Return type
None
Bonds
Description
The Bonds model describes pairwise potentials for 2 connected particles. Model creation occurs with a kwargs constructor as shown by equivalent operations below:
>>> bonds = mmelemental.models.forcefield.Bonds(
**{
"form": "Harmonic",
"params": {"spring": [2512.08, 2512.08]},
"lengths": [0.9572, 0.9572],
"connectivity": [[0,1,1.0], [0,2,1.0]],
}
)
This model can be used as an input to ForceField for an extended description of chemical bonds.
>>> ff = mmelemental.models.ForceField(
**{
"name": "water_ff",
"symbols": ["O", "H", "H"],
"charges": [-0.834, 0.417, 0.417],
"masses": [16.0, 1.008, 1.008],
"exclusions": "3",
"defs": ["OW", "HW", "HW"],
"nonbonded": nb,
"bonds": bonds,
}
)
>>> ff
ForceField(name='water_ff', form=['NonBonded', 'Bonds'], hash='236e292')
API
- class mmelemental.models.forcefield.bonded.bonds.Bonds(**data)[source]
- Parameters
form (str) – Name or form of the potential. See cls.supported_potentials for available potentials.
version (str, Optional) – Version of the force field this model stores. This field can be arbitrary.
params (Any) – Specific force field parameters model.
defs (List[str], Optional) – Particle definition. For atomic forcefields, this could be the atom type (e.g. HH31) or SMIRKS (OFF) representation.
extras (Dict[Any], Optional) – Additional information to bundle with the object. Use for schema development and scratch space.
lengths (Array) – Equilibrium bond lengths. Default unit is Angstroms.
lengths_units (str, Default: angstroms) – Equilibrium bond lengths unit.
connectivity (List[Tuple[Union[int, str], Union[int, str], float]], Optional) – Particle indices or names e.g. types for each bond and the bond order: (index1, index2, order).
- Return type
None
Angles
Description
The Angles model describes 3-body angle potentials for 3 connected particles. Model creation occurs with a kwargs constructor as shown by equivalent operations below:
>>> angles = mmelemental.models.forcefield.Angles(
**{
"form": "Harmonic",
"params": {"spring": [0.09565291598722435]},
"angles": [104.52],
"connectivity": [
[0, 1, 2],
],
}
)
This model can be used as an input to ForceField for an extended description of 3-angle terms in the force field.
>>> ff = mmelemental.models.ForceField(
**{
"name": "water_ff",
"symbols": ["O", "H", "H"],
"charges": [-0.834, 0.417, 0.417],
"masses": [16.0, 1.008, 1.008],
"exclusions": "3",
"defs": ["OW", "HW", "HW"],
"nonbonded": nb,
"bonds": bonds,
"angles": angles,
}
)
>>> ff
ForceField(name='water_ff', form=['NonBonded', 'Bonds', 'Angles'], hash='6bdf578')
API
- class mmelemental.models.forcefield.bonded.angles.Angles(**data)[source]
- Parameters
form (str) – Name or form of the potential. See cls.supported_potentials for available potentials.
version (str, Optional) – Version of the force field this model stores. This field can be arbitrary.
params (Any) – Specific force field parameters model.
defs (List[str], Optional) – Particle definition. For atomic forcefields, this could be the atom type (e.g. HH31) or SMIRKS (OFF) representation.
extras (Dict[Any], Optional) – Additional information to bundle with the object. Use for schema development and scratch space.
angles (Array) – Equilibrium angles. Default unit is degrees.
angles_units (str, Default: degrees) – Equilibrium angle units.
connectivity (List[Tuple[Union[int, str], Union[int, str], Union[int, str]]], Optional) – Particle indices for each angle.
- Return type
None
Dihedrals
Description
To be completed.
API
- class mmelemental.models.forcefield.bonded.dihedrals.Dihedrals(**data)[source]
- Parameters
form (str) – Name or form of the potential. See cls.supported_potentials for available potentials.
version (str, Optional) – Version of the force field this model stores. This field can be arbitrary.
params (Any) – Specific force field parameters model.
defs (List[str], Optional) – Particle definition. For atomic forcefields, this could be the atom type (e.g. HH31) or SMIRKS (OFF) representation.
extras (Dict[Any], Optional) – Additional information to bundle with the object. Use for schema development and scratch space.
angles (Array, Optional) – Equilibrium dihedral angles. Default unit is degrees.
angles_units (str, Default: degrees) – Equilibrium dihedral angle units.
connectivity (List[Tuple[Union[int, str], Union[int, str], Union[int, str], Union[int, str]]], Optional) – Particle indices for each dihedral angle.
weights (Array, Optional) – Something to consider later on? Ses CHARMM dihedral_style for LAMMPS.
- Return type
None
DihedralsImproper
Description
To be completed.
API
- class mmelemental.models.forcefield.bonded.dihedrals_improper.DihedralsImproper(**data)[source]
- Parameters
form (str) – Name or form of the potential. See cls.supported_potentials for available potentials.
version (str, Optional) – Version of the force field this model stores. This field can be arbitrary.
params (Any) – Specific force field parameters model.
defs (List[str], Optional) – Particle definition. For atomic forcefields, this could be the atom type (e.g. HH31) or SMIRKS (OFF) representation.
extras (Dict[Any], Optional) – Additional information to bundle with the object. Use for schema development and scratch space.
angles (Array, Optional) – Equilibrium dihedral angles. Default unit is degrees.
angles_units (str, Default: degrees) – Equilibrium improper dihedral angle units.
connectivity (List[Tuple[Union[int, str], Union[int, str], Union[int, str], Union[int, str]]], Optional) – Particle indices for each improper dihedral angle.
weights (Array, Optional) – Something to consider later on? Ses CHARMM dihedral_style for LAMMPS.
- Return type
None