OpenMM
|
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(NonbondedForce self) -> int More... | |
def | getNumExceptions |
getNumExceptions(NonbondedForce self) -> int More... | |
def | getNonbondedMethod |
getNonbondedMethod(NonbondedForce self) -> OpenMM::NonbondedForce::NonbondedMethod More... | |
def | setNonbondedMethod |
setNonbondedMethod(NonbondedForce self, OpenMM::NonbondedForce::NonbondedMethod method) More... | |
def | getCutoffDistance |
getCutoffDistance(NonbondedForce self) -> double More... | |
def | setCutoffDistance |
setCutoffDistance(NonbondedForce self, double distance) More... | |
def | getUseSwitchingFunction |
getUseSwitchingFunction(NonbondedForce self) -> bool More... | |
def | setUseSwitchingFunction |
setUseSwitchingFunction(NonbondedForce self, bool use) More... | |
def | getSwitchingDistance |
getSwitchingDistance(NonbondedForce self) -> double More... | |
def | setSwitchingDistance |
setSwitchingDistance(NonbondedForce self, double distance) More... | |
def | getReactionFieldDielectric |
getReactionFieldDielectric(NonbondedForce self) -> double More... | |
def | setReactionFieldDielectric |
setReactionFieldDielectric(NonbondedForce self, double dielectric) More... | |
def | getEwaldErrorTolerance |
getEwaldErrorTolerance(NonbondedForce self) -> double More... | |
def | setEwaldErrorTolerance |
setEwaldErrorTolerance(NonbondedForce self, double tol) More... | |
def | getPMEParameters |
getPMEParameters(NonbondedForce self) More... | |
def | setPMEParameters |
setPMEParameters(NonbondedForce self, double alpha, int nx, int ny, int nz) More... | |
def | addParticle |
addParticle(NonbondedForce self, double charge, double sigma, double epsilon) -> int More... | |
def | getParticleParameters |
getParticleParameters(NonbondedForce self, int index) More... | |
def | setParticleParameters |
setParticleParameters(NonbondedForce self, int index, double charge, double sigma, double epsilon) More... | |
def | addException |
addException(NonbondedForce self, int particle1, int particle2, double chargeProd, double sigma, double epsilon, bool replace=False) -> int addException(NonbondedForce self, int particle1, int particle2, double chargeProd, double sigma, double epsilon) -> int More... | |
def | getExceptionParameters |
getExceptionParameters(NonbondedForce self, int index) More... | |
def | setExceptionParameters |
setExceptionParameters(NonbondedForce self, int index, int particle1, int particle2, double chargeProd, double sigma, double epsilon) More... | |
def | createExceptionsFromBonds |
createExceptionsFromBonds(NonbondedForce self, vectorpairii bonds, double coulomb14Scale, double lj14Scale) More... | |
def | getUseDispersionCorrection |
getUseDispersionCorrection(NonbondedForce self) -> bool More... | |
def | setUseDispersionCorrection |
setUseDispersionCorrection(NonbondedForce self, bool useCorrection) More... | |
def | getReciprocalSpaceForceGroup |
getReciprocalSpaceForceGroup(NonbondedForce self) -> int More... | |
def | setReciprocalSpaceForceGroup |
setReciprocalSpaceForceGroup(NonbondedForce self, int group) More... | |
def | updateParametersInContext |
updateParametersInContext(NonbondedForce self, Context context) More... | |
def | addParticle_usingRVdw |
Add particle using elemetrary charge. More... | |
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. More... | |
def | __init__ |
init(OpenMM::NonbondedForce self) -> NonbondedForce init(OpenMM::NonbondedForce self, NonbondedForce other) -> NonbondedForce More... | |
def | __del__ |
del(OpenMM::NonbondedForce self) More... | |
Public Member Functions inherited from Force | |
def | __init__ |
def | __del__ |
del(OpenMM::Force self) More... | |
def | getForceGroup |
getForceGroup(Force self) -> int More... | |
def | setForceGroup |
setForceGroup(Force self, int group) More... | |
def | __copy__ |
def | __deepcopy__ |
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-Berthelot 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(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().
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.
When using a cutoff, by default Lennard-Jones interactions are sharply truncated at the cutoff distance. Optionally you can instead use a switching function to make the interaction smoothly go to zero over a finite distance range. To enable this, call setUseSwitchingFunction(). You must also call setSwitchingDistance() to specify the distance at which the interaction should begin to decrease. The switching distance must be less than the cutoff distance.
Another optional feature of this class (enabled by default) is to add a contribution to the energy which approximates the effect of all Lennard-Jones interactions beyond the cutoff in a periodic system. When running a simulation at constant pressure, this can improve the quality of the result. Call setUseDispersionCorrection() to set whether this should be used.
def __init__ | ( | self, | |
args | |||
) |
init(OpenMM::NonbondedForce self) -> NonbondedForce init(OpenMM::NonbondedForce self, NonbondedForce other) -> NonbondedForce
Create a NonbondedForce.
References simtk.openmm.openmm.stripUnits().
def __del__ | ( | self | ) |
del(OpenMM::NonbondedForce self)
References simtk.openmm.openmm.stripUnits().
def addException | ( | self, | |
args | |||
) |
addException(NonbondedForce self, int particle1, int particle2, double chargeProd, double sigma, double epsilon, bool replace=False) -> int addException(NonbondedForce 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. |
References simtk.openmm.openmm.stripUnits().
Referenced by NonbondedForce.addException_usingRMin().
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.
References NonbondedForce.addException().
Referenced by NonbondedForce.addParticle_usingRVdw().
def addParticle | ( | self, | |
args | |||
) |
addParticle(NonbondedForce 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-Berthelot 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 |
References simtk.openmm.openmm.stripUnits().
Referenced by NonbondedForce.addParticle_usingRVdw().
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
References NonbondedForce.addException_usingRMin(), AmoebaVdwForce.addParticle(), CustomGBForce.addParticle(), and NonbondedForce.addParticle().
def createExceptionsFromBonds | ( | self, | |
args | |||
) |
createExceptionsFromBonds(NonbondedForce 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 |
References simtk.openmm.openmm.stripUnits().
def getCutoffDistance | ( | self | ) |
getCutoffDistance(NonbondedForce 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.
References simtk.openmm.openmm.stripUnits().
def getEwaldErrorTolerance | ( | self | ) |
getEwaldErrorTolerance(NonbondedForce 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.
For PME calculations, if setPMEParameters() is used to set alpha to something other than 0, this value is ignored.
References simtk.openmm.openmm.stripUnits().
def getExceptionParameters | ( | self, | |
args | |||
) |
getExceptionParameters(NonbondedForce 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 |
References simtk.openmm.openmm.stripUnits().
def getNonbondedMethod | ( | self | ) |
getNonbondedMethod(NonbondedForce self) -> OpenMM::NonbondedForce::NonbondedMethod
Get the method used for handling long range nonbonded interactions.
References simtk.openmm.openmm.stripUnits().
def getNumExceptions | ( | self | ) |
getNumExceptions(NonbondedForce self) -> int
Get the number of special interactions that should be calculated differently from other interactions.
References simtk.openmm.openmm.stripUnits().
def getNumParticles | ( | self | ) |
getNumParticles(NonbondedForce self) -> int
Get the number of particles for which force field parameters have been defined.
References simtk.openmm.openmm.stripUnits().
def getParticleParameters | ( | self, | |
args | |||
) |
getParticleParameters(NonbondedForce 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 |
References simtk.openmm.openmm.stripUnits().
def getPMEParameters | ( | self | ) |
getPMEParameters(NonbondedForce self)
Get the parameters to use for PME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.
alpha | the separation parameter |
nx | the number of grid points along the X axis |
ny | the number of grid points along the Y axis |
nz | the number of grid points along the Z axis |
References simtk.openmm.openmm.stripUnits().
def getReactionFieldDielectric | ( | self | ) |
getReactionFieldDielectric(NonbondedForce self) -> double
Get the dielectric constant to use for the solvent in the reaction field approximation.
References simtk.openmm.openmm.stripUnits().
def getReciprocalSpaceForceGroup | ( | self | ) |
getReciprocalSpaceForceGroup(NonbondedForce 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.
References simtk.openmm.openmm.stripUnits().
def getSwitchingDistance | ( | self | ) |
getSwitchingDistance(NonbondedForce self) -> double
Get the distance at which the switching function begins to reduce the Lennard-Jones interaction. This must be less than the cutoff distance.
References simtk.openmm.openmm.stripUnits().
def getUseDispersionCorrection | ( | self | ) |
getUseDispersionCorrection(NonbondedForce 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.
References simtk.openmm.openmm.stripUnits().
def getUseSwitchingFunction | ( | self | ) |
getUseSwitchingFunction(NonbondedForce self) -> bool
Get whether a switching function is applied to the Lennard-Jones interaction. If the nonbonded method is set to NoCutoff, this option is ignored.
References simtk.openmm.openmm.stripUnits().
def setCutoffDistance | ( | self, | |
args | |||
) |
setCutoffDistance(NonbondedForce 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 |
References simtk.openmm.openmm.stripUnits().
def setEwaldErrorTolerance | ( | self, | |
args | |||
) |
setEwaldErrorTolerance(NonbondedForce self, double tol)
Set 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.
For PME calculations, if setPMEParameters() is used to set alpha to something other than 0, this value is ignored.
References simtk.openmm.openmm.stripUnits().
def setExceptionParameters | ( | self, | |
args | |||
) |
setExceptionParameters(NonbondedForce 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 |
References simtk.openmm.openmm.stripUnits().
def setNonbondedMethod | ( | self, | |
args | |||
) |
setNonbondedMethod(NonbondedForce self, OpenMM::NonbondedForce::NonbondedMethod method)
Set the method used for handling long range nonbonded interactions.
References simtk.openmm.openmm.stripUnits().
def setParticleParameters | ( | self, | |
args | |||
) |
setParticleParameters(NonbondedForce 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-Berthelot 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 |
References simtk.openmm.openmm.stripUnits().
def setPMEParameters | ( | self, | |
args | |||
) |
setPMEParameters(NonbondedForce self, double alpha, int nx, int ny, int nz)
Set the parameters to use for PME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.
alpha | the separation parameter |
nx | the number of grid points along the X axis |
ny | the number of grid points along the Y axis |
nz | the number of grid points along the Z axis |
References simtk.openmm.openmm.stripUnits().
def setReactionFieldDielectric | ( | self, | |
args | |||
) |
setReactionFieldDielectric(NonbondedForce self, double dielectric)
Set the dielectric constant to use for the solvent in the reaction field approximation.
References simtk.openmm.openmm.stripUnits().
def setReciprocalSpaceForceGroup | ( | self, | |
args | |||
) |
setReciprocalSpaceForceGroup(NonbondedForce 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. |
References simtk.openmm.openmm.stripUnits().
def setSwitchingDistance | ( | self, | |
args | |||
) |
setSwitchingDistance(NonbondedForce self, double distance)
Set the distance at which the switching function begins to reduce the Lennard-Jones interaction. This must be less than the cutoff distance.
References simtk.openmm.openmm.stripUnits().
def setUseDispersionCorrection | ( | self, | |
args | |||
) |
setUseDispersionCorrection(NonbondedForce 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.
References simtk.openmm.openmm.stripUnits().
def setUseSwitchingFunction | ( | self, | |
args | |||
) |
setUseSwitchingFunction(NonbondedForce self, bool use)
Set whether a switching function is applied to the Lennard-Jones interaction. If the nonbonded method is set to NoCutoff, this option is ignored.
References simtk.openmm.openmm.stripUnits().
def updateParametersInContext | ( | self, | |
args | |||
) |
updateParametersInContext(NonbondedForce self, Context context)
Update the particle and exception parameters in a Context to match those stored in this Force object. This method provides an efficient method to update certain parameters in an existing Context without needing to reinitialize it. Simply call setParticleParameters() and setExceptionParameters() to modify this object's parameters, then call updateParametersInState() to copy them over to the Context.
This method has several limitations. The only information it updates is the parameters of particles and exceptions. All other aspects of the Force (the nonbonded method, the cutoff distance, etc.) are unaffected and can only be changed by reinitializing the Context. Furthermore, only the chargeProd, sigma, and epsilon values of an exception can be changed; the pair of particles involved in the exception cannot change. Finally, this method cannot be used to add new particles or exceptions, only to change the parameters of existing ones.
References simtk.openmm.openmm.stripUnits().
|
static |
|
static |
|
static |
|
static |
|
static |
this |
Referenced by System.__init__().