OPLSPotentialBuilder¶
Included in QATK.Calculators.ForceField
- class OPLSPotentialBuilder(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, use_long_range_correction=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_untyped_molecules=None, add_missing_potentials=None)¶
This object builds a OPLS-AA potential using TremoloX.
- Parameters:
use_bonds – Controls if any bonding (stretching, angle, torsion, inversion and hydrogen bond) 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.
use_long_range_correction (bool) – Whether to use the long-range correction for the Lennard-Jones terms. Default: True
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_untyped_molecules (bool) – Whether whole bonded structures are ignored if they do not have potential atom types. Default: False
add_missing_potentials (bool) – Whether to add potentials for terms that are missing parameters. All types starting with the potential prefix (
UFF_,DREI_,OPLS_, etc.) will be considered as potential types. Default: False
- addCustomBondChargeIncrement(label_a, label_b, q_a, q_b)¶
Add a custom charge increments. Charge increments are used to calculate atomic charges. Atoms have their charge adjusted by the charge increments on the bonds attached to them. If opposite charges are added to each atom in the increment the whole configuration charge remains the same.
- Parameters:
label_a (str) – The first particle label. This can be either a type or a class.
label_b (str) – The second particle label. This can be either a type or a class.
q_a (PhysicalQuantity of type charge.) – The charge added to particle type_a.
q_b (PhysicalQuantity of type charge.) – The charge added to particle type_b.
- addCustomConstrainedAngle(class_a, class_b, class_c, theta_0)¶
Add a custom constrained angle between three bonded particle classes.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
theta_0 (PhysicalQuantity of type length.) – The equilibrium bond angle.
- addCustomConstrainedBond(class_a, class_b, r_0)¶
Add a custom constrained bond between two bonded particle classes.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
r_0 (PhysicalQuantity of type length.) – The equilibrium bond length.
- addCustomCosinePowerTorsion(class_a, class_b, class_c, class_d, k, c_list)¶
Add a custom power-series torsion potential between four bonded particle classes. The torsion potential is of the form:
\[V(\theta) = k \left[c_0 + \sum_{i=1} c_i \cos^i(\phi) \right]\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
class_d (str) – The fourth particle class.
k (PhysicalQuantity of type energy per square angle.) – The strength of the potential.
c_list (sequence of float) – The list of coefficients in the power series.
- addCustomCosineTorsion(class_a, class_b, class_c, class_d, k_list, n_list, delta_list)¶
Add a custom Cosine torsion between four bonded particle classes. The torsion potential is of the form:
\[V(\phi) = \sum_{i=1}^n k_i \left[1 + \cos(n_i \phi - \delta_i)\right]\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
class_d (str) – The fourth particle class.
k_list (Sequence of PhysicalQuantity of type energy.) – The size of the cosine potentials. This must have the same number of elements and n_list and delta_list.
n_list (Sequence of int.) – The period of each cosine function. This must have the same number of elements and k_list and delta_list.
delta_list (Sequence of PhysicalQuantity of type angle.) – The offset of the cosine potentials. This must have the same number of elements and k_list and n_list.
- addCustomDreidingHydrogenBond(donor_class, hydrogen_class, acceptor_class, a, b)¶
Add a custom Fourier improper torsion between four bonded particle types. The improper torsion potential is of the form:
\[V(\r, \theta) = a \left[5(r/b)^{12} - 6(r/b)^{10}\right] \cos^4 \theta\]The distance is measured between the hydrogen donor and acceptor, while the angle is measured with the hydrogen in the center. When the angle is less than 90 degrees the potential is zero.
- Parameters:
donor_class (str) – The donor particle class, which is bonded to the hydrogen.
hydrogen_class (str) – The hydrogen particle class.
acceptor_class (str) – The acceptor particle class, which is not bonded to the hydrogen.
a (PhysicalQuantity of type energy.) – The strength of the hydrogen bond potential.
b (PhysicalQuantity of type distance.) – The equilibrium hydrogen bond distance.
- addCustomFourierAngle(class_a, class_b, class_c, k, c_list)¶
Add a custom Fourier angle potential between three bonded particle classes. The angle potential is of the form:
\[V(\theta) = k \left[c_0 + \sum_{i=1}^4 c_i \cos(i \theta) \right]\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
k (PhysicalQuantity of type energy per square angle.) – The strength of the potential.
c_list (sequence of float) – The list of coefficients in the Fourier series.
- addCustomFourierImproperTorsion(class_a, class_b, class_c, class_d, k, c_list)¶
Add a custom Fourier improper torsion between four bonded particle classes. The improper torsion potential is of the form:
\[V(\chi) = \frac{k}{3} \sum_{i=0}^n \cos(i \chi)\]- Parameters:
class_a (str) – The central particle class.
class_b (str) – The first bonded particle class.
class_c (str) – The second bonded particle class.
class_d (str) – The third bonded particle class.
k (PhysicalQuantity of type energy per square angle.) – The strength of the potential.
c_list (sequence of float) – The list of coefficients in the Fourier series.
- addCustomFourierTorsion(class_a, class_b, class_c, class_d, k, c_list)¶
Add a custom Fourier torsion potential between four bonded particle classes. The torsion potential is of the form:
\[V(\theta) = k \left[c_0 + \sum_{i=1} c_i \cos(i \phi) \right]\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
class_d (str) – The fourth particle class.
k (PhysicalQuantity of type energy per square angle.) – The strength of the potential.
c_list (sequence of float) – The list of coefficients in the Fourier series.
- addCustomHarmonicAngle(class_a, class_b, class_c, k, theta_0)¶
Add a custom harmonic angle potential between three bonded particle classes. The angle potential is of the form:
\[V(\theta) = k(\theta - \theta_0)^2\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
k (PhysicalQuantity of type energy per square angle.) – The steepness of the harmonic potential.
theta_0 (PhysicalQuantity of type length.) – The equilibrium bond angle.
- addCustomHarmonicBond(class_a, class_b, k, r_0)¶
Add a custom harmonic bond potential between two bonded particle classes. The bond potential is of the form:
\[V(r) = k(r - r_0)^2\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
k (PhysicalQuantity of type energy per square distance.) – The steepness of the harmonic potential.
r_0 (PhysicalQuantity of type length.) – The equilibrium bond length.
- addCustomHarmonicCisTransTorsion(class_a, class_b, class_c, class_d, cis_k, trans_k, cis_torsion)¶
Add a custom harmonic torsion potential with optionally different parameters for torsions that are originally in the cis or trans orientation. This is most commonly applied to molecules that have double bonds. If the cis and trans parameters are the same, the same potential is used regardless of the initial torsion angle.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
class_d (str) – The fourth particle class.
cis_k (PhysicalQuantity of type energy per angle squared.) – The strength of the harmonic potential for a cis torsion.
trans_k (PhysicalQuantity of type energy per angle squared.) – The strength of the harmonic potential for a trans torsion.
cis_torsion (PhysicalQuantity of type angle.) – The equilibrium torsion angle of the cis configuration. The trans equilibrium angle is 180 degrees plus this angle.
- addCustomHarmonicImproperTorsion(class_a, class_b, class_c, class_d, k, psi_0)¶
Add a custom harmonic improper torsion potential between four bonded particle classes. The bond potential is of the form:
\[V(\psi) = k(\psi - \psi_0)^2\]- Parameters:
class_a (str) – The central particle class.
class_b (str) – The first bonded particle class.
class_c (str) – The second bonded particle class.
class_d (str) – The third bonded particle class.
k (PhysicalQuantity of type energy per square angle.) – The steepness of the harmonic potential.
psi_0 (PhysicalQuantity of type length.) – The equilibrium improper torsion angle.
- addCustomLennardJonesNonbond(type_a, type_b, epsilon, sigma)¶
Add a custom Lennard-Jones potential between two non-bonded atoms. These potentials are used instead of the potentials defined by the single atom parameters and combination rule.
- Parameters:
type_a (str) – The first particle type.
type_b (str) – The second particle type.
epsilon (PhysicalQuantity of type energy.) – The well depth of the potential.
sigma (PhysicalQuantity of type length.) – The distance where the potential is zero.
- addCustomOplsTorsion(class_a, class_b, class_c, class_d, k_1=PhysicalQuantity(0.0, eV), k_2=PhysicalQuantity(0.0, eV), k_3=PhysicalQuantity(0.0, eV), k_4=PhysicalQuantity(0.0, eV))¶
Add a custom OPLS torsion between four bonded particle classes. The torsion potential is of the form:
\[V(\phi) = \sum_{i=1}^4 \frac{k_i}{2}\left[1 + (-1)^{i+1} \cos(i \phi)\right]\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
class_d (str) – The fourth particle class.
k_1 (PhysicalQuantity of type energy.) – The first cosine potential size.
k_2 (PhysicalQuantity of type energy.) – The second cosine potential size.
k_3 (PhysicalQuantity of type energy.) – The third cosine potential size.
k_4 (PhysicalQuantity of type energy.) – The fourth cosine potential size.
- addCustomParticleType(type_label, class_label, description, particle_symbol, sigma, epsilon, charge, bond_number, type_filter=None, necessary_elements=None, definitions=None, overrides=None)¶
Add a new particle type to the potential.
- Parameters:
type_label (str) – The label for the new particle types. This is specified without the forcefield prefix.
class_label (str) – The class for the new particle type. Classes are used to group types together that use the same valence potentials such as bond or angle potentials.
description (str) – A description of the particle type.
particle_symbol (str) – The chemical symbol for the particle that the type refers to. Can be either an element symbol or a symbol from a UnitedAtom or Particle
sigma (PhysicalQuantity of type length.) – The Lennard-Jones sigma parameter for the type. This should correspond to where the Lennard-Jones potential is zero.
epsilon (PhysicalQuantity of type energy.) – The Lennard-Jones epsilon parameter for the type. This should correspond to the well depth in the potential.
charge (PhysicalQuantity of type charge.) – The charge of the particle type. This charge is assigned to the particle before charge increments are added to the particle charges.
bond_number (int) – The bond coordination number of the particle.
type_filter (str) – A label used to describe the chemical class that the type is found in. This can also be used to select types that generally work together.
necessary_elements (list of PeriodicTableElement) – A list of elements that must be present in the structure for the type to be found. This is used to speed atom typing by not searching for types for which the necessary atoms are not present. All types are assumed to require carbon, hydrogen and the element represented by the type.
definitions (str | list of str) – A SMARTS string defining the local bond graph of the atom type. This definition should be unique amongst the other user defined types. Multiple definitions, specified as a list of str, can also be given.
overrides (list of str) – A list of atom types that this type takes precedence over. This is used to resolve cases where an atom matches one ore more types, and allows some type definitions to be quite general. All user defined types are assumed to take precedence over pre-defined types.
- addCustomRepulsionShortRange(class_a, class_b, a)¶
Add a custom short-range repulsion potential between two bonded particle classes. The bond potential is of the form:
\[V(r) = \frac{a}{r^{12}}\]- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
a (PhysicalQuantity of type energy distance to the 12th power.) – The strength of the repulsion potential.
- addMissingPotentials()¶
- Returns:
Whether default potentials are added for terms that are missing explicit terms.
- Return type:
bool
- assignAtomCharges(configuration)¶
Assign the atomic partial charges to the given configuration. This only includes charges defined by the potential, which in some cases are zero.
- Parameters:
configuration (
MoleculeConfiguration|BulkConfiguration|DeviceConfiguration|SurfaceConfiguration) – The configuration to assign partial charges to.
- assignAtomTypes(configuration, search_types=None, overwrite=False, progress_widget=None)¶
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 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 to assign new types if the configuration already has a valid set of types. Default: False.
progress_widget (ForcefieldTypesWidget) – The progress bar / cancel button widget.
- checkAtomTypes(configuration)¶
Function that tests that each atom has only one valid potential type.
- Parameters:
configuration (AtomicConfiguration) – The configuration that is being tested.
- combinationRule()¶
- Returns:
The combination rule used for the Lennard-Jones potential.
- Return type:
str
- createCalculator(configuration, check_types=True, atomic_charges=None, charge_averaging=None)¶
Create and return a calculator with the potential.
- Parameters:
configuration (AtomicConfiguration) – The configuration the potential is to be used with.
check_types (bool) – Whether the types are checked for basic consistency. Default: True
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 (length)
- electrostaticScale()¶
- Returns:
Scaling factor for the electrostatic terms.
- Return type:
float
- electrostaticSmoothingLength()¶
- Returns:
Distance at which electrostatic terms are splined when using direct summation.
- Return type:
PhysicalQuantity (length)
- electrostaticSpmeCutoff()¶
- Returns:
The cutoff used in SPME electrostatic summation.
- Return type:
PhysicalQuantity (length)
- hBondCutoff()¶
- Returns:
The cutoff used in hydrogen bond summation.
- Return type:
PhysicalQuantity of type length
- hBondScale()¶
- Returns:
Scaling factor for the hydrogen bond terms.
- Return type:
float
- ignoreUntypedMolecules()¶
- Returns:
Whether bonded structures missing atom types are ignored in setting potentials.
- Return type:
bool
- includeBondAngles()¶
- Returns:
Whether bond angle terms are included in the potential set.
- Return type:
bool
- includeBondStretching()¶
- Returns:
Whether bond stretching terms are included in the potential set.
- Return type:
bool
- includeElectrostatic()¶
- Returns:
Whether electrostatic terms are included in the potential set.
- Return type:
bool
- includeHBonds()¶
- Returns:
Whether hydrogen bond terms are included in the potential set.
- Return type:
bool
- includeInversions()¶
- Returns:
Whether bond inversion terms are included in the potential set.
- Return type:
bool
- includeLennardJones()¶
- Returns:
Whether Lennard Jones terms are included in the potential set.
- Return type:
bool
- includeShortRange()¶
- Returns:
Whether short-range terms are included in the potential set.
- Return type:
bool
- includeTorsions()¶
- Returns:
Whether 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 the tag is defined in the potential.
- Return type:
bool
- lennardJonesCutoff()¶
- Returns:
Distance at which Lennard Jones terms are ignored.
- Return type:
PhysicalQuantity (length)
- lennardJonesSmoothingLength()¶
- Returns:
Distance at which Lennard Jones terms are splined.
- Return type:
PhysicalQuantity (length)
- oplsType(index, configuration)¶
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.
- Returns:
The type of the atom.
- Return type:
str
- potentialName()¶
- Returns:
The name used to identify the potential.
- 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.
- removeCustomAngle(class_a, class_b, class_c)¶
Remove a user-defined angle potential from the builder. When a custom potential overrides a standard potential, removing it restores the original standard potential parameters.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
- removeCustomBond(class_a, class_b)¶
Remove a user-defined bond potential from the builder. When a custom potential overrides a standard bond potential, removing it restores the original standard potential parameters.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
- removeCustomChargeIncrement(label_a, label_b)¶
Remove a user-defined charge increment from the builder. When a custom increment overrides a standard increment, removing it restores the original standard increment parameters.
- Parameters:
label_a – The first particle type.
type_b – The second particle type.
- removeCustomHydrogenBond(donor_class, hydrogen_class, acceptor_class)¶
Remove a user-defined hydrogen bond potential from the builder. When a custom potential overrides a standard potential, removing it restores the original standard potential parameters.
- Parameters:
donor_class (str) – The donor particle class, which is bonded to the hydrogen.
hydrogen_class (str) – The hydrogen particle class.
acceptor_class (str) – The acceptor particle class, which is not bonded to the hydrogen.
- removeCustomImproperTorsion(class_a, class_b, class_c, class_d)¶
Remove a user-defined improper torsion potential from the builder. When a custom potential overrides a standard potential, removing it restores the original standard potential parameters.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
class_d (str) – The fourth particle class.
- removeCustomNonbond(type_a, type_b)¶
Remove a user-defined Lennard-Jones potential from the builder. When a custom potential overrides a standard potential, removing it restores the original standard potential parameters.
- Parameters:
type_a (str) – The first particle type.
type_b (str) – The second particle type.
- removeCustomParticleType(type_label)¶
Remove a user-defined particle type from the builder. When a custom type overrides a standard particle type, removing it restores the original standard potential parameters.
- Parameters:
type_label (str) – The particle type.
- removeCustomShortRange(class_a, class_b)¶
Remove a user-defined short-range potential from the builder. When a custom potential overrides a standard potential, removing it restores the original standard potential parameters.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
- removeCustomTorsion(class_a, class_b, class_c, class_d)¶
Remove a user-defined torsion potential from the builder. In cases where the standard potential is over-written, this returns the standard potential.
- Parameters:
class_a (str) – The first particle class.
class_b (str) – The second particle class.
class_c (str) – The third particle class.
class_d (str) – The fourth particle class.
- scaleCoulomb()¶
- Returns:
The scaling factor for 1-4 interactions.
- Return type:
float
- scaleLennardJones()¶
- Returns:
The scaling factor for 1-4 interactions.
- Return type:
float
- scaleLennardJonesNonbond(type_a, type_b, epsilon_scale=1, sigma_scale=1)¶
Scale a Lennard-Jones potential by a given factor. This creates a user entry for this potential pair that over-writes the standard value.
- Parameters:
type_a (str) – The first particle type.
type_b (str) – The second particle type
epsilon_scale (float) – The scaling factor to apply to the well depth.
sigma_scale (float) – The scaling factor to apply to the particle size.
- tagPrefix()¶
- Returns:
The tag prefix used to identify the potential.
- Return type:
str
- totalCharge(configuration, verbose=True)¶
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.
- uniqueString()¶
Return a unique string representing the state of the object.
- useBonds()¶
- Returns:
Whether bond terms are included in the potential set.
- Return type:
bool
- useElectrostaticDirectSum()¶
- Returns:
Whether direct Coulomb sums are used for molecules.
- Return type:
bool
- useLongRangeCorrection()¶
- Returns:
Whether the long-range correction is used for the Lennard-Jones terms.
- Return type:
bool
Usage Examples¶
Build a bulk configuration of poly(vinyl chloride) plasticized with dimethyl sebacate using the OPLS potential.
# ---------------------------------------------
# Monomer And Included Molecules
# ---------------------------------------------
PVC_monomer = nlread('PVC_Monomer.hdf5')[-1]
dimthyl_sebacate = nlread('Dimethyl_Sebacate.hdf5')[-1]
# -------------------------------------------------------------
# OPLS Potential Definition
# -------------------------------------------------------------
opls_potential_builder = OPLSPotentialBuilder()
# Display the tags and charges placed on the dimethyl sebacate molecule
opls_potential_builder.printUsedTypes(dimthyl_sebacate)
opls_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,
)
bulk_configuration = polymer_builder.buildPolymerConfiguration()
# -------------------------------------------------------------
# Add OPLS Potential To Polymer Configuration
# -------------------------------------------------------------
calculator = opls_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_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-AA 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 overall charge of the molecule is zero.
Notes¶
The OPLSPotentialBuilder class enables the creation of a TremoloXCalculator that
implements the OPLS-AA 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 overall mathematical form of the OPLS-AA potential can be given as:
Here \(r\), \(\theta\), \(\phi\) and \(\chi\) represent inter-atomic and bond distances, angles, torsions and inversions respectively in the configuration.
To create a OPLS-AA calculator, an instance of the OPLSPotentialBuilder must first
be created. The constructor for this class has a number of options that controls the way in which
the calculator built. Individual types of interactions can be excluded by setting the include_
options to False By default all possible terms are included. All valence terms can also be
excluded by setting the argument use_bonds to False. This is appropriate for configurations
consisting of single atoms that are not bonded to each other.
As the OPLSPotentialBuilder class uses the bonding in the configuration to determine
appropriate potentials, bonds should be set on the configuration. This can be done with the
findBonds() method in either the MoleculeConfiguration or
BulkConfiguration class. This automatically assigns bonds between atoms that are
within the covalent radii of each other. Specific bonds can also be assigned with the
setBonds(), which takes a list of the indices of pairs of atoms that are
to be considered bonded.
There are then a number of options that control the non-bonding interactions. The parameter
lennard_jones_cutoff sets the distance beyond which Lennard-Jones interactions are truncated. To
avoid integration problems in molecular dynamics, Lennard-Jones interactions are brought to zero at
the cutoff using a spline function. The width of this spline can be set using the
lennard_jones_smoothing_length parameter. There are no specific hydrogen hydrogen bonding terms
defined within the OPLS-AA potential.
Electrostatic interactions are a little more complicated. By default the
OPLSPotentialBuilder adds a CoulombSPME electrostatic solver for
calculators to be used with a bulk configuration, and a CoulombDSF electrostatic solver
to calculators to be used with a molecule configuration. When adding an CoulombSPME
solver the cutoff used for calculating the real-space interactions can be set using the
electrostatic_spme_cutoff argument. Likewise the accuracy of the smoothed particle mesh Ewald
solver can be set using the electrostatic_accuracy argument. In the case where a
CoulombDSF solver is added to the potential, the cutoff used can be set using the
electrostatic_pairwise_cutoff argument. In cases where molecular electrostatics are to be
compared with values calculated in a bulk configuration, an explicit splined Coulomb potential
can be used for molecules instead of the damped shifted force potential. In this case a
CoulombN2Spline electrostatic solver can be added by using
use_electrostatic_direct_sum argument. As with the damped shifted force potential the cutoff is
controlled using the electrostatic_pairwise_cutoff argument, and like the Lennard-Jones
potential is brought to zero with a spline whose length is defined by the
electrostatic_smoothing_length argument. For all Coulomb solvers the total electrostatic
energies and forces can be scaled up or down using the electrostatic_scale argument.
Once the OPLSPotentialBuilder object is instantiated, a calculator is created using
the createCalculator method. This takes a configuration and returns the appropriate OPLS-AA
calculator for that configuration. To use electrostatic potentials in the calculator, atomic
partial charges have to be assigned to each atom. The OPLS-AA defines appropriate charges for each
atom type, which by default are added to the created calculator. It is also possible to add
pre-calculated charges from different sources using the atomic_charges argument. This argument
takes a list of charges that have the same ordering as the atoms in the configuration. The partial
atomic charges can also be averaged according to a few different schemes. This is controlled by the
charge_averaging parameter, which takes an NLFlag which indicates the type of averaging. The
flag AtomConnectivity averages charges on each atom based on the element and the other elements
it is bonded to. The flag AtomTypes averages according to the assigned OPLS-AA type. For the
calculator to recognize appropriate inversion terms, improper torsion indices also need to be set on
the configuration. By default this is done automatically, however existing improper torsions can be
preserved by setting the set_improper_torsion_indices to False. As the OPLS-AA potential
contains a number of different potentials, details of each potential used to create the calculator
can be printed out using the verbose argument, which is False by default.
When creating a OPLS-AA calculator for a configuration, the configuration must have the appropriate
tags for each atom type. Atom types in the OPLS-AA forcefield are grouped into classes. Atom types
in a class all have the same valence potential terms and individual Lennard-Jones and electrostatic
potentials. This reduces the number of parameters required to construct a transferable forcefield
while still allowing specific tailoring of the potentials to different molecules. Tags used by the
OPLS-AA forcefield all begin with the prefix OPLS_. The next section of tag gives the label of
the class of the atomic type, and the final section gives the number of the type in its class.
Appropriate OPLS tags can be assigned automatically using the method assignAtomTypes. This
method uses the bonding environment of each atom to estimate appropriate types. These tags are
then added to the input configuration. Type definitions have been added for atom types common to
polymers, although not all atom types have a definition. In cases where the
OPLSPotentialBuilder is unable to recognize the type of an atom, the tag OPLS_UNKNOWN
is added to that atom.
The OPLSPotentialBuilder contains a few methods to assist in manually assigning the correct types
to configurations in cases where the types cannot be automatically assigned. Atom types in the
OPLS-AA potential are generally defined for specific functional groups. As a result, atomic partial
charges for a given functional group usually add to zero or an integer charge for charged groups.
To discover the appropriate types for a configuration, the method printTotalTypes prints a
table of the available atomic types, with a short description of the kinds of atoms this type
represents. As they are defined together, types for different functional groups are often grouped
together. Types for just one class can also be printed using the type_class argument. The types
assigned on a specific configuration can also be displayed using the printUsedTypes method. As
functional groups generally have zero overall charge, correctly typed configurations should also
have zero overall charge. Calculating the total assigned charge to the configuration is therefore a
useful check to make sure the correct types have been added to the configuration. The total assigned
charge to a configuration using the OPLS forcefield can be calculated using the totalCharge
method.