This class implements nonbonded interactions between particles, including a Coulomb force to represent electrostatics and a Lennard-Jones force to represent van der Waals interactions. More...
Public Member Functions | |
def | getNumParticles |
getNumParticles(self) -> int | |
def | getNumExceptions |
getNumExceptions(self) -> int | |
def | getNonbondedMethod |
getNonbondedMethod(self) -> NonbondedMethod | |
def | setNonbondedMethod |
setNonbondedMethod(self, NonbondedMethod method) | |
def | getCutoffDistance |
getCutoffDistance(self) -> double | |
def | setCutoffDistance |
setCutoffDistance(self, double distance) | |
def | getReactionFieldDielectric |
getReactionFieldDielectric(self) -> double | |
def | setReactionFieldDielectric |
setReactionFieldDielectric(self, double dielectric) | |
def | getEwaldErrorTolerance |
getEwaldErrorTolerance(self) -> double | |
def | setEwaldErrorTolerance |
setEwaldErrorTolerance(self, double tol) | |
def | addParticle |
addParticle(self, double charge, double sigma, double epsilon) -> int | |
def | getParticleParameters |
getParticleParameters(self, int index) | |
def | setParticleParameters |
setParticleParameters(self, int index, double charge, double sigma, double epsilon) | |
def | addException |
addException(self, int particle1, int particle2, double chargeProd, double sigma, double epsilon, bool replace = False) -> int addException(self, int particle1, int particle2, double chargeProd, double sigma, double epsilon) -> int | |
def | getExceptionParameters |
getExceptionParameters(self, int index) | |
def | setExceptionParameters |
setExceptionParameters(self, int index, int particle1, int particle2, double chargeProd, double sigma, double epsilon) | |
def | createExceptionsFromBonds |
createExceptionsFromBonds(self, vectorpairii bonds, double coulomb14Scale, double lj14Scale) | |
def | getUseDispersionCorrection |
getUseDispersionCorrection(self) -> bool | |
def | setUseDispersionCorrection |
setUseDispersionCorrection(self, bool useCorrection) | |
def | getReciprocalSpaceForceGroup |
getReciprocalSpaceForceGroup(self) -> int | |
def | setReciprocalSpaceForceGroup |
setReciprocalSpaceForceGroup(self, int group) | |
def | addParticle_usingRVdw |
Add particle using elemetrary charge. | |
def | addException_usingRMin |
Add interaction exception using the product of the two atoms' elementary charges, rMin and epsilon, which is standard for AMBER force fields. | |
def | __init__ |
__init__(self) -> NonbondedForce __init__(self, NonbondedForce other) -> NonbondedForce | |
def | __del__ |
__del__(self) | |
Public Attributes | |
this | |
Static Public Attributes | |
NoCutoff = _openmm.NonbondedForce_NoCutoff | |
CutoffNonPeriodic = _openmm.NonbondedForce_CutoffNonPeriodic | |
CutoffPeriodic = _openmm.NonbondedForce_CutoffPeriodic | |
Ewald = _openmm.NonbondedForce_Ewald | |
PME = _openmm.NonbondedForce_PME |
This class implements nonbonded interactions between particles, including a Coulomb force to represent electrostatics and a Lennard-Jones force to represent van der Waals interactions.
It optionally supports periodic boundary conditions and cutoffs for long range interactions. Lennard-Jones interactions are calculated with the Lorentz-Bertelot combining rule: it uses the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles.
To use this class, create a NonbondedForce object, then call addParticle() once for each particle in the System to define its parameters. The number of particles for which you define nonbonded parameters must be exactly equal to the number of particles in the System, or else an exception will be thrown when you try to create a Context. After a particle has been added, you can modify its force field parameters by calling setParticleParameters().
NonbondedForce also lets you specify "exceptions", particular pairs of particles whose interactions should be computed based on different parameters than those defined for the individual particles. This can be used to completely exclude certain interactions from the force calculation, or to alter how they interact with each other.
Many molecular force fields omit Coulomb and Lennard-Jones interactions between particles separated by one or two bonds, while using modified parameters for those separated by three bonds (known as "1-4 interactions"). This class provides a convenience method for this case called createExceptionsFromBonds(). You pass to it a list of bonds and the scale factors to use for 1-4 interactions. It identifies all pairs of particles which are separated by 1, 2, or 3 bonds, then automatically creates exceptions for them.
def __del__ | ( | self | ) |
__del__(self)
Reimplemented from Force.
def __init__ | ( | self, | ||
args | ||||
) |
__init__(self) -> NonbondedForce __init__(self, NonbondedForce other) -> NonbondedForce
Create a NonbondedForce.
def addException | ( | self, | ||
args | ||||
) |
addException(self, int particle1, int particle2, double chargeProd, double sigma, double epsilon, bool replace = False) -> int addException(self, int particle1, int particle2, double chargeProd, double sigma, double epsilon) -> int
Add an interaction to the list of exceptions that should be calculated differently from other interactions. If chargeProd and epsilon are both equal to 0, this will cause the interaction to be completely omitted from force and energy calculations.
In many cases, you can use createExceptionsFromBonds() rather than adding each exception explicitly.
particle1 | the index of the first particle involved in the interaction | |
particle2 | the index of the second particle involved in the interaction | |
chargeProd | the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared | |
sigma | the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm | |
epsilon | the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol | |
replace | determines the behavior if there is already an exception for the same two particles. If true, the existing one is replaced. If false, an exception is thrown. |
def addException_usingRMin | ( | self, | ||
particle1, | ||||
particle2, | ||||
chargeProd, | ||||
rMin, | ||||
epsilon | ||||
) |
Add interaction exception using the product of the two atoms' elementary charges, rMin and epsilon, which is standard for AMBER force fields.
Note that rMin is the minimum energy point in the Lennard Jones potential. The conversion from sigma is: rMin = 2^1/6 * sigma.
def addParticle | ( | self, | ||
args | ||||
) |
addParticle(self, double charge, double sigma, double epsilon) -> int
Add the nonbonded force parameters for a particle. This should be called once for each particle in the System. When it is called for the i'th time, it specifies the parameters for the i'th particle. For calculating the Lennard-Jones interaction between two particles, the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles is used (the Lorentz-Bertelot combining rule).
charge | the charge of the particle, measured in units of the proton charge | |
sigma | the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm | |
epsilon | the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol |
def addParticle_usingRVdw | ( | self, | ||
charge, | ||||
rVDW, | ||||
epsilon | ||||
) |
Add particle using elemetrary charge.
Rvdw and epsilon, which is consistent with AMBER parameter file usage. Note that the sum of the radii of the two interacting atoms is the minimum energy point in the Lennard Jones potential and is often called rMin. The conversion from sigma follows: rVDW = 2^1/6 * sigma/2
def createExceptionsFromBonds | ( | self, | ||
args | ||||
) |
createExceptionsFromBonds(self, vectorpairii bonds, double coulomb14Scale, double lj14Scale)
Identify exceptions based on the molecular topology. Particles which are separated by one or two bonds are set to not interact at all, while pairs of particles separated by three bonds (known as "1-4 interactions") have their Coulomb and Lennard-Jones interactions reduced by a fixed factor.
bonds | the set of bonds based on which to construct exceptions. Each element specifies the indices of two particles that are bonded to each other. | |
coulomb14Scale | pairs of particles separated by three bonds will have the strength of their Coulomb interaction multiplied by this factor | |
lj14Scale | pairs of particles separated by three bonds will have the strength of their Lennard-Jones interaction multiplied by this factor |
def getCutoffDistance | ( | self | ) |
getCutoffDistance(self) -> double
Get the cutoff distance (in nm) being used for nonbonded interactions. If the NonbondedMethod in use is NoCutoff, this value will have no effect.
def getEwaldErrorTolerance | ( | self | ) |
getEwaldErrorTolerance(self) -> double
Get the error tolerance for Ewald summation. This corresponds to the fractional error in the forces which is acceptable. This value is used to select the reciprocal space cutoff and separation parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.
def getExceptionParameters | ( | self, | ||
args | ||||
) |
getExceptionParameters(self, int index)
Get the force field parameters for an interaction that should be calculated differently from others.
index | the index of the interaction for which to get parameters | |
particle1 | the index of the first particle involved in the interaction | |
particle2 | the index of the second particle involved in the interaction | |
chargeProd | the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared | |
sigma | the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm | |
epsilon | the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol |
def getNonbondedMethod | ( | self | ) |
getNonbondedMethod(self) -> NonbondedMethod
Get the method used for handling long range nonbonded interactions.
def getNumExceptions | ( | self | ) |
getNumExceptions(self) -> int
Get the number of special interactions that should be calculated differently from other interactions.
def getNumParticles | ( | self | ) |
getNumParticles(self) -> int
Get the number of particles for which force field parameters have been defined.
def getParticleParameters | ( | self, | ||
args | ||||
) |
getParticleParameters(self, int index)
Get the nonbonded force parameters for a particle.
index | the index of the particle for which to get parameters | |
charge | the charge of the particle, measured in units of the proton charge | |
sigma | the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm | |
epsilon | the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol |
def getReactionFieldDielectric | ( | self | ) |
getReactionFieldDielectric(self) -> double
Get the dielectric constant to use for the solvent in the reaction field approximation.
def getReciprocalSpaceForceGroup | ( | self | ) |
getReciprocalSpaceForceGroup(self) -> int
Get the force group that reciprocal space interactions for Ewald or PME are included in. This allows multiple time step integrators to evaluate direct and reciprocal space interactions at different intervals: getForceGroup() specifies the group for direct space, and getReciprocalSpaceForceGroup() specifies the group for reciprocal space. If this is -1 (the default value), the same force group is used for reciprocal space as for direct space.
def getUseDispersionCorrection | ( | self | ) |
getUseDispersionCorrection(self) -> bool
Get whether to add a contribution to the energy that approximately represents the effect of Lennard-Jones interactions beyond the cutoff distance. The energy depends on the volume of the periodic box, and is only applicable when periodic boundary conditions are used. When running simulations at constant pressure, adding this contribution can improve the quality of results.
def setCutoffDistance | ( | self, | ||
args | ||||
) |
setCutoffDistance(self, double distance)
Set the cutoff distance (in nm) being used for nonbonded interactions. If the NonbondedMethod in use is NoCutoff, this value will have no effect.
distance | the cutoff distance, measured in nm |
def setEwaldErrorTolerance | ( | self, | ||
args | ||||
) |
setEwaldErrorTolerance(self, double tol)
Get the error tolerance for Ewald summation. This corresponds to the fractional error in the forces which is acceptable. This value is used to select the reciprocal space cutoff and separation parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.
def setExceptionParameters | ( | self, | ||
args | ||||
) |
setExceptionParameters(self, int index, int particle1, int particle2, double chargeProd, double sigma, double epsilon)
Set the force field parameters for an interaction that should be calculated differently from others. If chargeProd and epsilon are both equal to 0, this will cause the interaction to be completely omitted from force and energy calculations.
index | the index of the interaction for which to get parameters | |
particle1 | the index of the first particle involved in the interaction | |
particle2 | the index of the second particle involved in the interaction | |
chargeProd | the scaled product of the atomic charges (i.e. the strength of the Coulomb interaction), measured in units of the proton charge squared | |
sigma | the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm | |
epsilon | the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol |
def setNonbondedMethod | ( | self, | ||
args | ||||
) |
setNonbondedMethod(self, NonbondedMethod method)
Set the method used for handling long range nonbonded interactions.
def setParticleParameters | ( | self, | ||
args | ||||
) |
setParticleParameters(self, int index, double charge, double sigma, double epsilon)
Set the nonbonded force parameters for a particle. When calculating the Lennard-Jones interaction between two particles, it uses the arithmetic mean of the sigmas and the geometric mean of the epsilons for the two interacting particles (the Lorentz-Bertelot combining rule).
index | the index of the particle for which to set parameters | |
charge | the charge of the particle, measured in units of the proton charge | |
sigma | the sigma parameter of the Lennard-Jones potential (corresponding to the van der Waals radius of the particle), measured in nm | |
epsilon | the epsilon parameter of the Lennard-Jones potential (corresponding to the well depth of the van der Waals interaction), measured in kJ/mol |
def setReactionFieldDielectric | ( | self, | ||
args | ||||
) |
setReactionFieldDielectric(self, double dielectric)
Set the dielectric constant to use for the solvent in the reaction field approximation.
def setReciprocalSpaceForceGroup | ( | self, | ||
args | ||||
) |
setReciprocalSpaceForceGroup(self, int group)
Set the force group that reciprocal space interactions for Ewald or PME are included in. This allows multiple time step integrators to evaluate direct and reciprocal space interactions at different intervals: setForceGroup() specifies the group for direct space, and setReciprocalSpaceForceGroup() specifies the group for reciprocal space. If this is -1 (the default value), the same force group is used for reciprocal space as for direct space.
group | the group index. Legal values are between 0 and 31 (inclusive), or -1 to use the same force group that is specified for direct space. |
def setUseDispersionCorrection | ( | self, | ||
args | ||||
) |
setUseDispersionCorrection(self, bool useCorrection)
Set whether to add a contribution to the energy that approximately represents the effect of Lennard-Jones interactions beyond the cutoff distance. The energy depends on the volume of the periodic box, and is only applicable when periodic boundary conditions are used. When running simulations at constant pressure, adding this contribution can improve the quality of results.
CutoffNonPeriodic = _openmm.NonbondedForce_CutoffNonPeriodic [static] |
CutoffPeriodic = _openmm.NonbondedForce_CutoffPeriodic [static] |
Ewald = _openmm.NonbondedForce_Ewald [static] |
NoCutoff = _openmm.NonbondedForce_NoCutoff [static] |
PME = _openmm.NonbondedForce_PME [static] |