OPLSMinPotentialBuilder¶
- class OPLSMinPotentialBuilder(use_bonds=None, include_bond_stretching=None, include_bond_angles=None, include_torsions=None, include_inversions=None, include_lennard_jones=None, include_electrostatic=None, lennard_jones_cutoff=None, lennard_jones_smoothing_length=None, electrostatic_pairwise_cutoff=None, electrostatic_smoothing_length=None, electrostatic_spme_cutoff=None, electrostatic_accuracy=None, electrostatic_scale=None, use_electrostatic_direct_sum=None, ignore_missing_types=None)¶
This object builds a OPLS-Min potential using TremoloX.
- Parameters:
use_bonds (bool) – Controls if any bonding (stretching, angle, torsion and inversion) terms are to be included. Default: True.
include_bond_stretching (bool) – Controls if bonding stretching terms are included. Default: use_bonds.
include_bond_angles (bool) – Controls if angle bending terms are included. Default: use_bonds.
include_torsions (bool) – Controls if bond torsion terms are included. Default: use_bonds.
include_inversions (bool) – Controls if bond inversion torsion terms are included. Default: use_bonds.
include_lennard_jones (bool) – Controls if Lennard-Jones terms are included. Default: use_bonds.
include_electrostatic (bool) – Controls if electrostatic terms are included. Default: use_bonds.
lennard_jones_cutoff (PhysicalQuantity of type distance) – Distance beyond which the Lennard-Jones terms are truncated. Default: 10.0 * Angstrom.
lennard_jones_smoothing_length (PhysicalQuantity of type distance) – Distance over which the Lennard-Jones terms are brought to zero at the cutoff. Default: 2.0 * Angstrom.
electrostatic_pairwise_cutoff (PhysicalQuantity of type distance) – Distance at which the pairwise electrostatic terms are brought to zero in pairwise summation. Default: 12 * Angstrom
electrostatic_smoothing_length (PhysicalQuantity of type distance) – Distance over which the electrostatic terms are brought to zero at the cutoff when using the direct electrostatic sum. Default: 2.0 * Angstrom
electrostatic_spme_cutoff (PhysicalQuantity of type distance) – Distance at which the pairwise electrostatic terms are brought to zero in SPME summation. Default: 7.5 * Angstrom
electrostatic_accuracy (float) – Sets the relative accuracy of the SPME summation. Default: 0.0001
electrostatic_scale (float) – Set the relative strength of the electrostatic terms. Default: 1.0
use_electrostatic_direct_sum (bool) – Use the direct Coulomb sum to calculate the electrostatic terms in MoleculeConfigurations. Default: False
ignore_missing_types (bool) – Whether or not molecules missing OPLS types are ignored. Default: False
- assignAtomTypes(configuration, search_types=None, overwrite=False)¶
Automatically assign OPLS atom tags to a configuration.
- Parameters:
configuration (AtomicConfiguration) – The configuration on which to assign tags.
search_types (list of str | None) – List of types to assign onto the configuration. If specified only these types will be assigned, with OPLS_UNKNOWN assigned to atoms that don’t match any of the given types. The default argument None means that all available types are assigned. Default: None.
overwrite (bool) – Whether or not to assign new types if the configuration already has a valid set of types. Default: False
- checkAtomTypes(configuration, add_missing_terms=False)¶
Function that tests that each atom has only one valid OPLS type.
- Parameters:
configuration (AtomicConfiguration) – The configuration that is being tested.
add_missing_terms (bool) – Whether or not to add terms for missing OPLS types. Default: False
- createCalculator(configuration, check_types=True, verbose=False, atomic_charges=None, charge_averaging=None)¶
Create and return a calculator with the OPLS-Min potential.
- Parameters:
configuration (AtomicConfiguration) – The configuration the potential is to be used with.
check_types (bool) – Flag as to whether or not the correctness of the types are checked. Default: True
verbose (bool) – Whether or not to print the terms assigned in the potential set. Default: False
atomic_charges (PhysicalQuantity of type charge) – Custom atomic charges to be used in the potential.
charge_averaging (NLFlag | None) – The type of averaging of charges to be used, if any. Accepted values are None, AtomConnectivity to average over atoms bonded to the same elements and AtomTypes to average over atoms of the same forcefield type. Default: None
- Returns:
The calculator containing the potential.
- Return type:
- electrostaticAccuracy()¶
- Returns:
The relative accuracy of the electrostatic summation for SPME.
- Return type:
float
- electrostaticPairwiseCutoff()¶
- Returns:
The cutoff used in pairwise electrostatic summation.
- Return type:
PhysicalQuantity of type length
- electrostaticScale()¶
- Returns:
Scaling factor for the electrostatic terms.
- Return type:
float
- electrostaticSmoothingLength()¶
- Returns:
Distance over which electrostatic terms are splined when using direct summation.
- Return type:
PhysicalQuantity of type length
- electrostaticSpmeCutoff()¶
- Returns:
The cutoff used in SPME electrostatic summation.
- Return type:
PhysicalQuantity of type length
- ignoreMissingTypes()¶
- Returns:
Whether or not missing types are ignored on molecules.
- Return type:
bool
- includeBondAngles()¶
- Returns:
Whether or not bond angle terms are included in the potential set.
- Return type:
bool
- includeBondStretching()¶
- Returns:
Whether or not bond stretching terms are included in the potential set.
- Return type:
bool
- includeElectrostatic()¶
- Returns:
Whether or not electrostatic terms are included in the potential set.
- Return type:
bool
- includeInversions()¶
- Returns:
Whether or not bond inversion terms are included in the potential set.
- Return type:
bool
- includeLennardJones()¶
- Returns:
Whether or not Lennard Jones terms are included in the potential set.
- Return type:
bool
- includeTorsions()¶
- Returns:
Whether or not bond torsion terms are included in the potential set.
- Return type:
bool
- isOplsType(type_label)¶
Function that tests if a type is in the potential.
- Parameters:
type_label (str) – The type being tested.
- Returns:
Whether or not the tag is defined in the potential.
- Return type:
bool
- lennardJonesCutoff()¶
- Returns:
Distance beyond which Lennard Jones terms are truncated.
- Return type:
PhysicalQuantity of type length
- lennardJonesSmoothingLength()¶
- Returns:
Distance over which Lennard Jones terms are splined.
- Return type:
PhysicalQuantity of type length
- nlinfo()¶
- Returns:
The nlinfo.
- Return type:
dict
- oplsType(index, configuration, add_missing_terms=False)¶
Gives the OPLS type on an atom in a configuration, ignoring other tags.
- Parameters:
index (int) – Index of the atom for which the type is sought.
configuration (AtomicConfiguration) – The configuration containing the atom.
add_missing_terms (bool) – Whether or not to add terms for missing OPLS types. Default: False
- Returns:
The type of the atom.
- Return type:
str
- printTotalTypes(type_class=None)¶
Print a table of available atom types. This provides a list of types and descriptions which can be useful when doing manual type assignment.
- Parameters:
type_class (str | None) – Display only atom types with this class label. None displays all types. Default: None
- printUsedTypes(configuration)¶
Print a table of the types used on the configuration.
- Parameters:
configuration (AtomicConfiguration) – The configuration whose types are being displayed.
- totalCharge(configuration, verbose=False)¶
Return the total charge of the configuration according to the atom types on the molecule.
- Parameters:
configuration (AtomicConfiguration) – The configuration whose charges are being reported.
verbose (bool) – Flag as to whether or not additional information about charges is printed. Default: False
- Returns:
The total charge of the configuration according to the atom types on the molecule.
- Return type:
PhysicalQuantity with units of charge.
- typePartialCharge(type_label)¶
Return the atomic partial charge associated with the given OPLS type label.
- Parameters:
type_label (str) – The type label for the OPLS atom type.
- Returns:
The atom charge associated with the type label.
- Return type:
PhysicalQuantity of type charge
- uniqueString()¶
Return a unique string representing the state of the object.
- useBonds()¶
- Returns:
Whether or not bond terms are included in the potential set.
- Return type:
bool
- useElectrostaticDirectSum()¶
- Returns:
Whether or not direct Coulomb sums are used for molecules.
- Return type:
bool
Usage Examples¶
Build a bulk configuration of poly(vinyl chloride) plasticized with dimethyl sebacate using the OPLS-Min potential.
# ---------------------------------------------
# Monomer And Included Molecules
# ---------------------------------------------
PVC_monomer = nlread('PVC_Monomer.hdf5')[-1]
dimthyl_sebacate = nlread('Dimethyl_Sebacate.hdf5')[-1]
# -------------------------------------------------------------
# OPLS Potential Definition
# -------------------------------------------------------------
opls_min_potential_builder = OPLSMinPotentialBuilder()
# Display the tags and charges placed on the dimethyl sebacate molecule
opls_min_potential_builder.printUsedTypes(dimthyl_sebacate)
opls_min_potential_builder.totalCharge(dimthyl_sebacate, verbose=True)
# ---------------------------------------------
# Polymer Sequence
# ---------------------------------------------
polymer_sequence = PolymerSequence(
monomer_configurations=[PVC_monomer],
number_of_monomers=20,
number_of_chains=20,
tactic_ratio=0.5,
)
# -------------------------------------------------------------
# Polymer Monte Carlo Builder
# -------------------------------------------------------------
polymer_builder = PolymerMonteCarloBuilder(
polymer_sequence=[polymer_sequence],
density=1.4000*gram/cm**3,
included_molecules=[(dimthyl_sebacate, 4)],
monte_carlo_temperature=300.0*Kelvin,
angle_sampling_points=20,
repack_molecules=False,
potential_builder=opls_min_potential_builder
)
bulk_configuration = polymer_builder.buildPolymerConfiguration()
# -------------------------------------------------------------
# Add OPLS Potential To Polymer Configuration
# -------------------------------------------------------------
calculator = opls_min_potential_builder.createCalculator(bulk_configuration)
bulk_configuration.setCalculator(calculator)
nlsave('PVC_DMS.hdf5', bulk_configuration)
# -------------------------------------------------------------
# Force Capped Equilibration
# -------------------------------------------------------------
equilibration_method = ForceCappedEquilibration(
temperature=300.00*Kelvin,
md_steps_per_force_capped_simulation=40000,
md_time_step=0.50*fs,
force_capped_simulations=4,
starting_factor=1.040,
ending_factor=0.800,
fixed_indices=None,
)
bulk_configuration = equilibration_method.runEquilibration(bulk_configuration)
nlsave('PVC_DMS.hdf5', bulk_configuration)
OPLS-Min_PVC_DMS_Example.py
PVC_Monomer.hdf5
Dimethyl_Sebacate.hdf5
The given example above demonstrates both how to create a calculator that uses the OPLS-Min potential and also how to assign tags for the dimethyl sebacate molecule that is included in the polymer melt. Both the atom types and charges are printed. This shows both the appropriate types that have been selected and also that the default atom types don’t have charge.
Notes¶
The OPLSMinPotentialBuilder class enables the creation of a TremoloXCalculator
that
implements the OPLS potential of Mayo, Olafson and Goddard[1]. This is a
bonded valence forcefield that represents the potential energy of a molecule as a collection of
simple functions based on bond lengths, bond angles, torsion angles, inversion angles and
inter-atomic distances.
The OPLSMinPotentialBuilder class is very similar to the OPLSPotentialBuilder
class.
Both implement different versions of the OPLS potential. In the OPLSPotentialBuilder
atom types have assigned atomic partial charges. This results in a large number of atom types
required to capture many molecular charge polarizations. In the OPLSMinPotentialBuilder
no charges are associated with each atom type, allowing each atom class to only have one type.
The significantly simplifies assigning the correct atom types, and is most useful in calculations
where atomic partial charges different to the ones in the OPLSPotentialBuilder
are
required.
Atomic partial charges for the configuration can be added when creating the calculator, or
estimated using charge equilibration. Both the OPLSMinPotentialBuilder and
OPLSPotentialBuilder
use the same potential parameters for the bond lengths, angles,
torsions and improper torsions in the configuration. There are however small differences in some
Lennard-Jones non-bonding parameters, resulting from only having one atom type per class. The
tagging scheme used by the OPLSMinPotentialBuilder broadly follows that used in the
OPLSPotentialBuilder
. Atom tags start with the prefix OPLS_
, then contain a string
representation of the atom class, and then a number representing the atom type number. In the
OPLSMinPotentialBuilder there is only one atom type per class, so this number is 0
in all
atom types.