CustomHbondForce Class Reference

This class supports a wide variety of energy functions used to represent hydrogen bonding. More...

Inheritance diagram for CustomHbondForce:
Force

List of all members.

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 getNumFunctions
 getNumFunctions(self) -> int
def getEnergyFunction
 getEnergyFunction(self) -> string
def setEnergyFunction
 setEnergyFunction(self, string energy)
def getNonbondedMethod
 getNonbondedMethod(self) -> NonbondedMethod
def setNonbondedMethod
 setNonbondedMethod(self, NonbondedMethod method)
def getCutoffDistance
 getCutoffDistance(self) -> double
def setCutoffDistance
 setCutoffDistance(self, double distance)
def addPerDonorParameter
 addPerDonorParameter(self, string name) -> int
def getPerDonorParameterName
 getPerDonorParameterName(self, int index) -> string
def setPerDonorParameterName
 setPerDonorParameterName(self, int index, string name)
def addPerAcceptorParameter
 addPerAcceptorParameter(self, string name) -> int
def getPerAcceptorParameterName
 getPerAcceptorParameterName(self, int index) -> string
def setPerAcceptorParameterName
 setPerAcceptorParameterName(self, int index, string name)
def addGlobalParameter
 addGlobalParameter(self, string name, double defaultValue) -> int
def getGlobalParameterName
 getGlobalParameterName(self, int index) -> string
def setGlobalParameterName
 setGlobalParameterName(self, int index, string name)
def getGlobalParameterDefaultValue
 getGlobalParameterDefaultValue(self, int index) -> double
def setGlobalParameterDefaultValue
 setGlobalParameterDefaultValue(self, int index, double defaultValue)
def addDonor
 addDonor(self, int d1, int d2, int d3, vectord parameters) -> int
def getDonorParameters
 getDonorParameters(self, int index)
def setDonorParameters
 setDonorParameters(self, int index, int d1, int d2, int d3, vectord parameters)
def addAcceptor
 addAcceptor(self, int a1, int a2, int a3, vectord parameters) -> int
def getAcceptorParameters
 getAcceptorParameters(self, int index)
def setAcceptorParameters
 setAcceptorParameters(self, int index, int a1, int a2, int a3, vectord parameters)
def addExclusion
 addExclusion(self, int donor, int acceptor) -> int
def getExclusionParticles
 getExclusionParticles(self, int index)
def setExclusionParticles
 setExclusionParticles(self, int index, int donor, int acceptor)
def addFunction
 addFunction(self, string name, vectord values, double min, double max) -> int
def getFunctionParameters
 getFunctionParameters(self, int index)
def setFunctionParameters
 setFunctionParameters(self, int index, string name, vectord values, double min, double max)
def __init__
 __init__(self, string energy) -> CustomHbondForce __init__(self, CustomHbondForce other) -> CustomHbondForce
def __del__
 __del__(self)

Public Attributes

 this

Static Public Attributes

 NoCutoff = _openmm.CustomHbondForce_NoCutoff
 CutoffNonPeriodic = _openmm.CustomHbondForce_CutoffNonPeriodic
 CutoffPeriodic = _openmm.CustomHbondForce_CutoffPeriodic

Detailed Description

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().

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, step. All trigonometric functions are defined in radians, and log is the natural logarithm. step(x) = 0 if x is less than 0, 1 otherwise.

In addition, you can call addFunction() to define a new function based on tabulated values. You specify a vector of values, and a natural spline is created from them. That function can then appear in the expression.


Member Function Documentation

def __del__ (   self  ) 

__del__(self)

Reimplemented from Force.

def __init__ (   self,
  args 
)

__init__(self, string energy) -> CustomHbondForce __init__(self, CustomHbondForce other) -> CustomHbondForce

Create a CustomHbondForce.

Parameters:
energy an algebraic expression giving the interaction energy between a donor 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, int a1, int a2, int a3, vectord parameters) -> int

Add an acceptor group to the force

Parameters:
a1 the index of the first particle for this acceptor group
a2 the index of the second particle for this acceptor group. If the group only includes one particle, this must be -1.
a3 the index of the third particle for this acceptor group. If the group includes less than three particles, this must be -1.
parameters the list of per-acceptor parameter values for the new acceptor
def addDonor (   self,
  args 
)

addDonor(self, int d1, int d2, int d3, vectord parameters) -> int

Add a donor group to the force

Parameters:
d1 the index of the first particle for this donor group
d2 the index of the second particle for this donor group. If the group only includes one particle, this must be -1.
d3 the index of the third particle for this donor group. If the group includes less than three particles, this must be -1.
parameters the list of per-donor parameter values for the new donor
def addExclusion (   self,
  args 
)

addExclusion(self, int donor, int acceptor) -> int

Add a donor-acceptor pair to the list of interactions that should be excluded.

Parameters:
donor the index of the donor to exclude
acceptor the index of the acceptor to exclude
def addFunction (   self,
  args 
)

addFunction(self, string name, vectord values, double min, double max) -> int

Add a tabulated function that may appear in the energy expression.

Parameters:
name the name of the function as it appears in expressions
values the tabulated values of the function f(x) at uniformly spaced values of x between min and max. The function is assumed to be zero for x < min or x > max.
min the value of the independent variable corresponding to the first element of values
max the value of the independent variable corresponding to the last element of values
def addGlobalParameter (   self,
  args 
)

addGlobalParameter(self, string name, double defaultValue) -> int

Add a new global parameter that the interaction may depend on.

Parameters:
name the name of the parameter
defaultValue the default value of the parameter
def addPerAcceptorParameter (   self,
  args 
)

addPerAcceptorParameter(self, string name) -> int

Add a new per-acceptor parameter that the interaction may depend on.

Parameters:
name the name of the parameter
def addPerDonorParameter (   self,
  args 
)

addPerDonorParameter(self, string name) -> int

Add a new per-donor parameter that the interaction may depend on.

Parameters:
name the name of the parameter
def getAcceptorParameters (   self,
  args 
)

getAcceptorParameters(self, int index)

Get the properties of an acceptor group.

Parameters:
index the index of the acceptor group to get
a1 the index of the first particle for this acceptor group
a2 the index of the second particle for this acceptor group. If the group only includes one particle, this will be -1.
a3 the index of the third particle for this acceptor group. If the group includes less than three particles, this will be -1.
parameters the list of per-acceptor parameter values for the new acceptor
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,
  args 
)

getDonorParameters(self, int index)

Get the properties of a donor group.

Parameters:
index the index of the donor group to get
d1 the index of the first particle for this donor group
d2 the index of the second particle for this donor group. If the group only includes one particle, this will be -1.
d3 the index of the third particle for this donor group. If the group includes less than three particles, this will be -1.
parameters the list of per-donor parameter values for the new donor
def getEnergyFunction (   self  ) 

getEnergyFunction(self) -> string

Get the algebraic expression that gives the interaction energy between a donor and an acceptor

def getExclusionParticles (   self,
  args 
)

getExclusionParticles(self, int index)

Get the donor and acceptor in a pair whose interaction should be excluded.

Parameters:
index the index of the exclusion for which to get donor and acceptor indices
particle1 the index of the donor
particle2 the index of the acceptor
def getFunctionParameters (   self,
  args 
)

getFunctionParameters(self, int index)

Get the parameters for a tabulated function that may appear in the energy expression.

Parameters:
index the index of the function for which to get parameters
name the name of the function as it appears in expressions
values the tabulated values of the function f(x) at uniformly spaced values of x between min and max. The function is assumed to be zero for x < min or x > max.
min the value of the independent variable corresponding to the first element of values
max the value of the independent variable corresponding to the last element of values
def getGlobalParameterDefaultValue (   self,
  args 
)

getGlobalParameterDefaultValue(self, int index) -> double

Get the default value of a global parameter.

Parameters:
index the index of the parameter for which to get the default value
def getGlobalParameterName (   self,
  args 
)

getGlobalParameterName(self, int index) -> string

Get the name of a global parameter.

Parameters:
index the index of the parameter for which to get the name
def getNonbondedMethod (   self  ) 

getNonbondedMethod(self) -> 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 getPerAcceptorParameterName (   self,
  args 
)

getPerAcceptorParameterName(self, int index) -> string

Get the name of a per-acceptor parameter.

Parameters:
index the index of the parameter for which to get the name
def getPerDonorParameterName (   self,
  args 
)

getPerDonorParameterName(self, int index) -> string

Get the name of a per-donor parameter.

Parameters:
index the index of the parameter for which to get the name
def setAcceptorParameters (   self,
  args 
)

setAcceptorParameters(self, int index, int a1, int a2, int a3, vectord parameters)

Set the properties of an acceptor group.

Parameters:
index the index of the acceptor group to set
a1 the index of the first particle for this acceptor group
a2 the index of the second particle for this acceptor group. If the group only includes one particle, this must be -1.
a3 the index of the third particle for this acceptor group. If the group includes less than three particles, this must be -1.
parameters the list of per-acceptor parameter values for the new acceptor
def setCutoffDistance (   self,
  args 
)

setCutoffDistance(self, double 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.

Parameters:
distance the cutoff distance, measured in nm
def setDonorParameters (   self,
  args 
)

setDonorParameters(self, int index, int d1, int d2, int d3, vectord parameters)

Set the properties of a donor group.

Parameters:
index the index of the donor group to set
d1 the index of the first particle for this donor group
d2 the index of the second particle for this donor group. If the group only includes one particle, this must be -1.
d3 the index of the third particle for this donor group. If the group includes less than three particles, this must be -1.
parameters the list of per-donor parameter values for the new donor
def setEnergyFunction (   self,
  args 
)

setEnergyFunction(self, string energy)

Set the algebraic expression that gives the interaction energy between a donor and an acceptor

def setExclusionParticles (   self,
  args 
)

setExclusionParticles(self, int index, int donor, int acceptor)

Get the donor and acceptor in a pair whose interaction should be excluded.

Parameters:
index the index of the exclusion for which to get donor and acceptor indices
particle1 the index of the donor
particle2 the index of the acceptor
def setFunctionParameters (   self,
  args 
)

setFunctionParameters(self, int index, string name, vectord values, double min, double max)

Set the parameters for a tabulated function that may appear in algebraic expressions.

Parameters:
index the index of the function for which to set parameters
name the name of the function as it appears in expressions
values the tabulated values of the function f(x) at uniformly spaced values of x between min and max. The function is assumed to be zero for x < min or x > max.
min the value of the independent variable corresponding to the first element of values
max the value of the independent variable corresponding to the last element of values
def setGlobalParameterDefaultValue (   self,
  args 
)

setGlobalParameterDefaultValue(self, int index, double defaultValue)

Set the default value of a global parameter.

Parameters:
index the index of the parameter for which to set the default value
name the default value of the parameter
def setGlobalParameterName (   self,
  args 
)

setGlobalParameterName(self, int index, string name)

Set the name of a global parameter.

Parameters:
index the index of the parameter for which to set the name
name the name of the parameter
def setNonbondedMethod (   self,
  args 
)

setNonbondedMethod(self, NonbondedMethod method)

Set the method used for handling long range nonbonded interactions.

def setPerAcceptorParameterName (   self,
  args 
)

setPerAcceptorParameterName(self, int index, string name)

Set the name of a per-acceptor parameter.

Parameters:
index the index of the parameter for which to set the name
name the name of the parameter
def setPerDonorParameterName (   self,
  args 
)

setPerDonorParameterName(self, int index, string name)

Set the name of a per-donor parameter.

Parameters:
index the index of the parameter for which to set the name
name the name of the parameter

Member Data Documentation

CutoffNonPeriodic = _openmm.CustomHbondForce_CutoffNonPeriodic [static]
CutoffPeriodic = _openmm.CustomHbondForce_CutoffPeriodic [static]
NoCutoff = _openmm.CustomHbondForce_NoCutoff [static]

The documentation for this class was generated from the following file:

Generated by  doxygen 1.6.2