OpenMM
|
This class supports a wide variety of energy functions used to represent hydrogen bonding. More...
Public Member Functions | |
def | getNumDonors |
getNumDonors(self) -> int | |
def | getNumAcceptors |
getNumAcceptors(self) -> int | |
def | getNumExclusions |
getNumExclusions(self) -> int | |
def | getNumPerDonorParameters |
getNumPerDonorParameters(self) -> int | |
def | getNumPerAcceptorParameters |
getNumPerAcceptorParameters(self) -> int | |
def | getNumGlobalParameters |
getNumGlobalParameters(self) -> int | |
def | getNumTabulatedFunctions |
getNumTabulatedFunctions(self) -> int | |
def | getNumFunctions |
getNumFunctions(self) -> int | |
def | getEnergyFunction |
getEnergyFunction(self) -> std::string const & | |
def | setEnergyFunction |
Set the algebraic expression that gives the interaction energy between a donor and an acceptor. | |
def | getNonbondedMethod |
getNonbondedMethod(self) -> OpenMM::CustomHbondForce::NonbondedMethod | |
def | setNonbondedMethod |
Set the method used for handling long range nonbonded interactions. | |
def | getCutoffDistance |
getCutoffDistance(self) -> double | |
def | setCutoffDistance |
Set the cutoff distance (in nm) being used. | |
def | addPerDonorParameter |
addPerDonorParameter(self, name) -> int | |
def | getPerDonorParameterName |
getPerDonorParameterName(self, index) -> std::string const & | |
def | setPerDonorParameterName |
Set the name of a per-donor parameter. | |
def | addPerAcceptorParameter |
addPerAcceptorParameter(self, name) -> int | |
def | getPerAcceptorParameterName |
getPerAcceptorParameterName(self, index) -> std::string const & | |
def | setPerAcceptorParameterName |
Set the name of a per-acceptor parameter. | |
def | addGlobalParameter |
addGlobalParameter(self, name, defaultValue) -> int | |
def | getGlobalParameterName |
getGlobalParameterName(self, index) -> std::string const & | |
def | setGlobalParameterName |
Set the name of a global parameter. | |
def | getGlobalParameterDefaultValue |
getGlobalParameterDefaultValue(self, index) -> double | |
def | setGlobalParameterDefaultValue |
Set the default value of a global parameter. | |
def | addDonor |
addDonor(self, d1, d2, d3, parameters) -> int addDonor(self, d1, d2, d3) -> int | |
def | getDonorParameters |
Get the properties of a donor group. | |
def | setDonorParameters |
Set the properties of a donor group. | |
def | addAcceptor |
addAcceptor(self, a1, a2, a3, parameters) -> int addAcceptor(self, a1, a2, a3) -> int | |
def | getAcceptorParameters |
Get the properties of an acceptor group. | |
def | setAcceptorParameters |
Set the properties of an acceptor group. | |
def | addExclusion |
addExclusion(self, donor, acceptor) -> int | |
def | getExclusionParticles |
Get the donor and acceptor in a pair whose interaction should be excluded. | |
def | setExclusionParticles |
Get the donor and acceptor in a pair whose interaction should be excluded. | |
def | addTabulatedFunction |
addTabulatedFunction(self, name, function) -> int | |
def | getTabulatedFunction |
getTabulatedFunction(self, index) -> TabulatedFunction getTabulatedFunction(self, index) -> TabulatedFunction | |
def | getTabulatedFunctionName |
getTabulatedFunctionName(self, index) -> std::string const & | |
def | addFunction |
addFunction(self, name, values, min, max) -> int | |
def | getFunctionParameters |
Get the parameters for a tabulated function that may appear in the energy expression. | |
def | setFunctionParameters |
Set the parameters for a tabulated function that may appear in the energy expression. | |
def | updateParametersInContext |
Update the per-donor and per-acceptor parameters in a Context to match those stored in this Force object. | |
def | usesPeriodicBoundaryConditions |
usesPeriodicBoundaryConditions(self) -> bool | |
def | __init__ |
__init__(self, energy) -> CustomHbondForce __init__(self, other) -> CustomHbondForce | |
Public Attributes | |
this | |
Static Public Attributes | |
NoCutoff = _openmm.CustomHbondForce_NoCutoff | |
CutoffNonPeriodic = _openmm.CustomHbondForce_CutoffNonPeriodic | |
CutoffPeriodic = _openmm.CustomHbondForce_CutoffPeriodic |
This class supports a wide variety of energy functions used to represent hydrogen bonding.
It computes interactions between "donor" particle groups and "acceptor" particle groups, where each group may include up to three particles. Typically a donor group consists of a hydrogen atom and the atoms it is bonded to, and an acceptor group consists of a negatively charged atom and the atoms it is bonded to.
We refer to the particles in a donor group as d1, d2 and d3, and the particles in an acceptor group as a1, a2, and a3. For each donor and each acceptor, CustomHbondForce evaluates a user supplied algebraic expression to determine the interaction energy. The expression may depend on arbitrary distances, angles, and dihedral angles defined by any of the six particles involved. The function distance(p1, p2) is the distance between the particles p1 and p2 (where "p1" and "p2" should be replaced by the names of the actual particles to calculate the distance between), angle(p1, p2, p3) is the angle formed by the three specified particles, and dihedral(p1, p2, p3, p4) is the dihedral angle formed by the four specified particles.
The expression also may involve tabulated functions, and may depend on arbitrary global, per-donor, and per-acceptor parameters. It also optionally supports periodic boundary conditions and cutoffs for long range interactions.
To use this class, create a CustomHbondForce object, passing an algebraic expression to the constructor that defines the interaction energy between each donor and acceptor. Then call addPerDonorParameter() to define per-donor parameters, addPerAcceptorParameter() to define per-acceptor parameters, and addGlobalParameter() to define global parameters. The values of per-donor and per-acceptor parameters are specified as part of the system definition, while values of global parameters may be modified during a simulation by calling Context::setParameter().
Next, call addDonor() and addAcceptor() to define donors and acceptors and specify their parameter values. After a donor or acceptor has been added, you can modify its parameters by calling setDonorParameters() or setAcceptorParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().
CustomHbondForce also lets you specify "exclusions", particular combinations of donors and acceptors whose interactions should be omitted from force and energy calculations. This is most often used for particles that are bonded to each other.
As an example, the following code creates a CustomHbondForce that implements a simple harmonic potential to keep the distance between a1 and d1, and the angle formed by a1-d1-d2, near ideal values:
CustomHbondForce* force = new CustomHbondForce("k*(distance(a1,d1)-r0)^2*(angle(a1,d1,d2)-theta0)^2");
This force depends on three parameters: k, r0, and theta0. The following code defines these as per-donor parameters:
force->addPerDonorParameter("k"); force->addPerDonorParameter("r0"); force->addPerDonorParameter("theta0");
Expressions may involve the operators + (add), - (subtract), * (multiply), / (divide), and ^ (power), and the following functions: sqrt, exp, log, sin, cos, sec, csc, tan, cot, asin, acos, atan, sinh, cosh, tanh, erf, erfc, min, max, abs, floor, ceil, step, delta, select. All trigonometric functions are defined in radians, and log is the natural logarithm. step(x) = 0 if x is less than 0, 1 otherwise. delta(x) = 1 if x is 0, 0 otherwise. select(x,y,z) = z if x = 0, y otherwise.
In addition, you can call addTabulatedFunction() to define a new function based on tabulated values. You specify the function by creating a TabulatedFunction object. That function can then appear in the expression.
def __init__ | ( | self, | |
args | |||
) |
__init__(self, energy) -> CustomHbondForce __init__(self, other) -> CustomHbondForce
Create a CustomHbondForce.
energy | (string) an algebraic expression giving the interaction energy between a donor and an acceptor as a function of inter-particle distances, angles, and dihedrals, as well as any global, per-donor, and per-acceptor parameters |
def addAcceptor | ( | self, | |
args | |||
) |
addAcceptor(self, a1, a2, a3, parameters) -> int addAcceptor(self, a1, a2, a3) -> int
Add an acceptor group to the force
a1 | (int) the index of the first particle for this acceptor group |
a2 | (int) the index of the second particle for this acceptor group. If the group only includes one particle, this must be -1. |
a3 | (int) the index of the third particle for this acceptor group. If the group includes less than three particles, this must be -1. |
parameters | (vector< double >) the list of per-acceptor parameter values for the new acceptor |
def addDonor | ( | self, | |
args | |||
) |
addDonor(self, d1, d2, d3, parameters) -> int addDonor(self, d1, d2, d3) -> int
Add a donor group to the force
d1 | (int) the index of the first particle for this donor group |
d2 | (int) the index of the second particle for this donor group. If the group only includes one particle, this must be -1. |
d3 | (int) the index of the third particle for this donor group. If the group includes less than three particles, this must be -1. |
parameters | (vector< double >) the list of per-donor parameter values for the new donor |
def addExclusion | ( | self, | |
donor, | |||
acceptor | |||
) |
addExclusion(self, donor, acceptor) -> int
Add a donor-acceptor pair to the list of interactions that should be excluded.
donor | (int) the index of the donor to exclude |
acceptor | (int) the index of the acceptor to exclude |
def addFunction | ( | self, | |
name, | |||
values, | |||
min, | |||
max | |||
) |
addFunction(self, name, values, min, max) -> int
Add a tabulated function that may appear in the energy expression.
def addGlobalParameter | ( | self, | |
name, | |||
defaultValue | |||
) |
addGlobalParameter(self, name, defaultValue) -> int
Add a new global parameter that the interaction may depend on.
name | (string) the name of the parameter |
defaultValue | (double) the default value of the parameter |
def addPerAcceptorParameter | ( | self, | |
name | |||
) |
addPerAcceptorParameter(self, name) -> int
Add a new per-acceptor parameter that the interaction may depend on.
name | (string) the name of the parameter |
def addPerDonorParameter | ( | self, | |
name | |||
) |
addPerDonorParameter(self, name) -> int
Add a new per-donor parameter that the interaction may depend on.
name | (string) the name of the parameter |
def addTabulatedFunction | ( | self, | |
name, | |||
function | |||
) |
addTabulatedFunction(self, name, function) -> int
Add a tabulated function that may appear in the energy expression.
name | (string) the name of the function as it appears in expressions |
function | (TabulatedFunction *) a TabulatedFunction object defining the function. The TabulatedFunction should have been created on the heap with the "new" operator. The Force takes over ownership of it, and deletes it when the Force itself is deleted. |
def getAcceptorParameters | ( | self, | |
index | |||
) |
Get the properties of an acceptor group.
index | (int) the index of the acceptor group to get |
def getCutoffDistance | ( | self | ) |
getCutoffDistance(self) -> double
Get the cutoff distance (in nm) being used. All interactions for which the distance between d1 and a1 is greater than the cutoff will be ignored. If the NonbondedMethod in use is NoCutoff, this value will have no effect.
def getDonorParameters | ( | self, | |
index | |||
) |
Get the properties of a donor group.
index | (int) the index of the donor group to get |
def getEnergyFunction | ( | self | ) |
getEnergyFunction(self) -> std::string const &
Get the algebraic expression that gives the interaction energy between a donor and an acceptor
def getExclusionParticles | ( | self, | |
index | |||
) |
Get the donor and acceptor in a pair whose interaction should be excluded.
index | (int) the index of the exclusion for which to get donor and acceptor indices |
def getFunctionParameters | ( | self, | |
index | |||
) |
Get the parameters for a tabulated function that may appear in the energy expression.
def getGlobalParameterDefaultValue | ( | self, | |
index | |||
) |
getGlobalParameterDefaultValue(self, index) -> double
Get the default value of a global parameter.
index | (int) the index of the parameter for which to get the default value |
def getGlobalParameterName | ( | self, | |
index | |||
) |
getGlobalParameterName(self, index) -> std::string const &
Get the name of a global parameter.
index | (int) the index of the parameter for which to get the name |
def getNonbondedMethod | ( | self | ) |
getNonbondedMethod(self) -> OpenMM::CustomHbondForce::NonbondedMethod
Get the method used for handling long range nonbonded interactions.
def getNumAcceptors | ( | self | ) |
getNumAcceptors(self) -> int
Get the number of acceptors for which force field parameters have been defined.
def getNumDonors | ( | self | ) |
getNumDonors(self) -> int
Get the number of donors for which force field parameters have been defined.
def getNumExclusions | ( | self | ) |
getNumExclusions(self) -> int
Get the number of donor-acceptor pairs whose interactions should be excluded.
def getNumFunctions | ( | self | ) |
getNumFunctions(self) -> int
Get the number of tabulated functions that have been defined.
def getNumGlobalParameters | ( | self | ) |
getNumGlobalParameters(self) -> int
Get the number of global parameters that the interaction depends on.
def getNumPerAcceptorParameters | ( | self | ) |
getNumPerAcceptorParameters(self) -> int
Get the number of per-acceptor parameters that the interaction depends on.
def getNumPerDonorParameters | ( | self | ) |
getNumPerDonorParameters(self) -> int
Get the number of per-donor parameters that the interaction depends on.
def getNumTabulatedFunctions | ( | self | ) |
getNumTabulatedFunctions(self) -> int
Get the number of tabulated functions that have been defined.
def getPerAcceptorParameterName | ( | self, | |
index | |||
) |
getPerAcceptorParameterName(self, index) -> std::string const &
Get the name of a per-acceptor parameter.
index | (int) the index of the parameter for which to get the name |
def getPerDonorParameterName | ( | self, | |
index | |||
) |
getPerDonorParameterName(self, index) -> std::string const &
Get the name of a per-donor parameter.
index | (int) the index of the parameter for which to get the name |
def getTabulatedFunction | ( | self, | |
args | |||
) |
getTabulatedFunction(self, index) -> TabulatedFunction getTabulatedFunction(self, index) -> TabulatedFunction
Get a reference to a tabulated function that may appear in the energy expression.
index | (int) the index of the function to get |
def getTabulatedFunctionName | ( | self, | |
index | |||
) |
getTabulatedFunctionName(self, index) -> std::string const &
Get the name of a tabulated function that may appear in the energy expression.
index | (int) the index of the function to get |
def setAcceptorParameters | ( | self, | |
args | |||
) |
Set the properties of an acceptor group.
index | (int) the index of the acceptor group to set |
a1 | (int) the index of the first particle for this acceptor group |
a2 | (int) the index of the second particle for this acceptor group. If the group only includes one particle, this must be -1. |
a3 | (int) the index of the third particle for this acceptor group. If the group includes less than three particles, this must be -1. |
parameters | (vector< double >) the list of per-acceptor parameter values for the acceptor |
def setCutoffDistance | ( | self, | |
distance | |||
) |
Set the cutoff distance (in nm) being used.
All interactions for which the distance between d1 and a1 is greater than the cutoff will be ignored. If the NonbondedMethod in use is NoCutoff, this value will have no effect.
distance | (double) the cutoff distance, measured in nm |
def setDonorParameters | ( | self, | |
args | |||
) |
Set the properties of a donor group.
index | (int) the index of the donor group to set |
d1 | (int) the index of the first particle for this donor group |
d2 | (int) the index of the second particle for this donor group. If the group only includes one particle, this must be -1. |
d3 | (int) the index of the third particle for this donor group. If the group includes less than three particles, this must be -1. |
parameters | (vector< double >) the list of per-donor parameter values for the donor |
def setEnergyFunction | ( | self, | |
energy | |||
) |
Set the algebraic expression that gives the interaction energy between a donor and an acceptor.
def setExclusionParticles | ( | self, | |
index, | |||
donor, | |||
acceptor | |||
) |
Get the donor and acceptor in a pair whose interaction should be excluded.
index | (int) the index of the exclusion for which to get donor and acceptor indices |
donor | (int) the index of the donor |
acceptor | (int) the index of the acceptor |
def setFunctionParameters | ( | self, | |
index, | |||
name, | |||
values, | |||
min, | |||
max | |||
) |
Set the parameters for a tabulated function that may appear in the energy expression.
def setGlobalParameterDefaultValue | ( | self, | |
index, | |||
defaultValue | |||
) |
Set the default value of a global parameter.
index | (int) the index of the parameter for which to set the default value |
defaultValue | (double) the default value of the parameter |
def setGlobalParameterName | ( | self, | |
index, | |||
name | |||
) |
Set the name of a global parameter.
index | (int) the index of the parameter for which to set the name |
name | (string) the name of the parameter |
def setNonbondedMethod | ( | self, | |
method | |||
) |
Set the method used for handling long range nonbonded interactions.
def setPerAcceptorParameterName | ( | self, | |
index, | |||
name | |||
) |
Set the name of a per-acceptor parameter.
index | (int) the index of the parameter for which to set the name |
name | (string) the name of the parameter |
def setPerDonorParameterName | ( | self, | |
index, | |||
name | |||
) |
Set the name of a per-donor parameter.
index | (int) the index of the parameter for which to set the name |
name | (string) the name of the parameter |
def updateParametersInContext | ( | self, | |
context | |||
) |
Update the per-donor and per-acceptor 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 setDonorParameters() and setAcceptorParameters() to modify this object's parameters, then call updateParametersInContext() to copy them over to the Context.
This method has several limitations. The only information it updates is the values of per-donor and per-acceptor parameters. All other aspects of the Force (the energy function, nonbonded method, cutoff distance, etc.) are unaffected and can only be changed by reinitializing the Context. The set of particles involved in a donor or acceptor cannot be changed, nor can new donors or acceptors be added.
def usesPeriodicBoundaryConditions | ( | self | ) |
usesPeriodicBoundaryConditions(self) -> bool
Returns whether or not this force makes use of periodic boundary conditions.
Reimplemented from Force.
CutoffNonPeriodic = _openmm.CustomHbondForce_CutoffNonPeriodic [static] |
CutoffPeriodic = _openmm.CustomHbondForce_CutoffPeriodic [static] |
NoCutoff = _openmm.CustomHbondForce_NoCutoff [static] |